1 /* 2 Objects to manage the interactions between the mesh data structures and the algebraic objects 3 */ 4 #if !defined(__PETSCDM_H) 5 #define __PETSCDM_H 6 #include <petscmat.h> 7 #include <petscdmtypes.h> 8 #include <petscfetypes.h> 9 10 PETSC_EXTERN PetscErrorCode DMInitializePackage(void); 11 12 PETSC_EXTERN PetscClassId DM_CLASSID; 13 14 /*J 15 DMType - String with the name of a PETSc DM 16 17 Level: beginner 18 19 .seealso: DMSetType(), DM 20 J*/ 21 typedef const char* DMType; 22 #define DMDA "da" 23 #define DMCOMPOSITE "composite" 24 #define DMSLICED "sliced" 25 #define DMSHELL "shell" 26 #define DMPLEX "plex" 27 #define DMCARTESIAN "cartesian" 28 #define DMREDUNDANT "redundant" 29 #define DMPATCH "patch" 30 #define DMMOAB "moab" 31 #define DMCIRCUIT "circuit" 32 33 PETSC_EXTERN const char *const DMBoundaryTypes[]; 34 PETSC_EXTERN PetscFunctionList DMList; 35 PETSC_EXTERN PetscBool DMRegisterAllCalled; 36 PETSC_EXTERN PetscErrorCode DMCreate(MPI_Comm,DM*); 37 PETSC_EXTERN PetscErrorCode DMClone(DM,DM*); 38 PETSC_EXTERN PetscErrorCode DMSetType(DM, DMType); 39 PETSC_EXTERN PetscErrorCode DMGetType(DM, DMType *); 40 PETSC_EXTERN PetscErrorCode DMRegister(const char[],PetscErrorCode (*)(DM)); 41 PETSC_EXTERN PetscErrorCode DMRegisterAll(void); 42 PETSC_EXTERN PetscErrorCode DMRegisterDestroy(void); 43 44 PETSC_EXTERN PetscErrorCode DMView(DM,PetscViewer); 45 PETSC_EXTERN PetscErrorCode DMLoad(DM,PetscViewer); 46 PETSC_EXTERN PetscErrorCode DMDestroy(DM*); 47 PETSC_EXTERN PetscErrorCode DMCreateGlobalVector(DM,Vec*); 48 PETSC_EXTERN PetscErrorCode DMCreateLocalVector(DM,Vec*); 49 PETSC_EXTERN PetscErrorCode DMGetLocalVector(DM,Vec *); 50 PETSC_EXTERN PetscErrorCode DMRestoreLocalVector(DM,Vec *); 51 PETSC_EXTERN PetscErrorCode DMGetGlobalVector(DM,Vec *); 52 PETSC_EXTERN PetscErrorCode DMRestoreGlobalVector(DM,Vec *); 53 PETSC_EXTERN PetscErrorCode DMClearGlobalVectors(DM); 54 PETSC_EXTERN PetscErrorCode DMGetNamedGlobalVector(DM,const char*,Vec*); 55 PETSC_EXTERN PetscErrorCode DMRestoreNamedGlobalVector(DM,const char*,Vec*); 56 PETSC_EXTERN PetscErrorCode DMGetNamedLocalVector(DM,const char*,Vec*); 57 PETSC_EXTERN PetscErrorCode DMRestoreNamedLocalVector(DM,const char*,Vec*); 58 PETSC_EXTERN PetscErrorCode DMGetLocalToGlobalMapping(DM,ISLocalToGlobalMapping*); 59 PETSC_EXTERN PetscErrorCode DMCreateFieldIS(DM,PetscInt*,char***,IS**); 60 PETSC_EXTERN PetscErrorCode DMGetBlockSize(DM,PetscInt*); 61 PETSC_EXTERN PetscErrorCode DMCreateColoring(DM,ISColoringType,ISColoring*); 62 PETSC_EXTERN PetscErrorCode DMCreateMatrix(DM,Mat*); 63 PETSC_EXTERN PetscErrorCode DMSetMatrixPreallocateOnly(DM,PetscBool); 64 PETSC_EXTERN PetscErrorCode DMCreateInterpolation(DM,DM,Mat*,Vec*); 65 PETSC_EXTERN PetscErrorCode DMCreateInjection(DM,DM,VecScatter*); 66 PETSC_EXTERN PetscErrorCode DMGetWorkArray(DM,PetscInt,PetscDataType,void*); 67 PETSC_EXTERN PetscErrorCode DMRestoreWorkArray(DM,PetscInt,PetscDataType,void*); 68 PETSC_EXTERN PetscErrorCode DMRefine(DM,MPI_Comm,DM*); 69 PETSC_EXTERN PetscErrorCode DMCoarsen(DM,MPI_Comm,DM*); 70 PETSC_EXTERN PetscErrorCode DMRefineHierarchy(DM,PetscInt,DM[]); 71 PETSC_EXTERN PetscErrorCode DMCoarsenHierarchy(DM,PetscInt,DM[]); 72 PETSC_EXTERN PetscErrorCode DMCoarsenHookAdd(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,Mat,Vec,Mat,DM,void*),void*); 73 PETSC_EXTERN PetscErrorCode DMRefineHookAdd(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,Mat,DM,void*),void*); 74 PETSC_EXTERN PetscErrorCode DMRestrict(DM,Mat,Vec,Mat,DM); 75 PETSC_EXTERN PetscErrorCode DMInterpolate(DM,Mat,DM); 76 PETSC_EXTERN PetscErrorCode DMSetFromOptions(DM); 77 PETSC_STATIC_INLINE PetscErrorCode DMViewFromOptions(DM A,const char prefix[],const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,prefix,name);} 78 79 PETSC_EXTERN PetscErrorCode DMSetUp(DM); 80 PETSC_EXTERN PetscErrorCode DMCreateInterpolationScale(DM,DM,Mat,Vec*); 81 PETSC_EXTERN PetscErrorCode DMCreateAggregates(DM,DM,Mat*); 82 PETSC_EXTERN PetscErrorCode DMGlobalToLocalHookAdd(DM,PetscErrorCode (*)(DM,Vec,InsertMode,Vec,void*),PetscErrorCode (*)(DM,Vec,InsertMode,Vec,void*),void*); 83 PETSC_EXTERN PetscErrorCode DMGlobalToLocalBegin(DM,Vec,InsertMode,Vec); 84 PETSC_EXTERN PetscErrorCode DMGlobalToLocalEnd(DM,Vec,InsertMode,Vec); 85 PETSC_EXTERN PetscErrorCode DMLocalToGlobalBegin(DM,Vec,InsertMode,Vec); 86 PETSC_EXTERN PetscErrorCode DMLocalToGlobalEnd(DM,Vec,InsertMode,Vec); 87 PETSC_EXTERN PetscErrorCode DMLocalToLocalBegin(DM,Vec,InsertMode,Vec); 88 PETSC_EXTERN PetscErrorCode DMLocalToLocalEnd(DM,Vec,InsertMode,Vec); 89 PETSC_EXTERN PetscErrorCode DMConvert(DM,DMType,DM*); 90 91 PETSC_EXTERN PetscErrorCode DMGetCoordinateDM(DM,DM*); 92 PETSC_EXTERN PetscErrorCode DMSetCoordinateDM(DM,DM); 93 PETSC_EXTERN PetscErrorCode DMGetCoordinateSection(DM,PetscSection*); 94 PETSC_EXTERN PetscErrorCode DMSetCoordinateSection(DM,PetscSection); 95 PETSC_EXTERN PetscErrorCode DMGetCoordinates(DM,Vec*); 96 PETSC_EXTERN PetscErrorCode DMSetCoordinates(DM,Vec); 97 PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocal(DM,Vec*); 98 PETSC_EXTERN PetscErrorCode DMSetCoordinatesLocal(DM,Vec); 99 PETSC_EXTERN PetscErrorCode DMLocatePoints(DM,Vec,IS*); 100 PETSC_EXTERN PetscErrorCode DMGetPeriodicity(DM,const PetscReal**,const PetscReal**); 101 PETSC_EXTERN PetscErrorCode DMSetPeriodicity(DM,const PetscReal[],const PetscReal[]); 102 103 /* block hook interface */ 104 PETSC_EXTERN PetscErrorCode DMSubDomainHookAdd(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,VecScatter,VecScatter,DM,void*),void*); 105 PETSC_EXTERN PetscErrorCode DMSubDomainRestrict(DM,VecScatter,VecScatter,DM); 106 107 PETSC_EXTERN PetscErrorCode DMSetOptionsPrefix(DM,const char []); 108 PETSC_EXTERN PetscErrorCode DMSetVecType(DM,VecType); 109 PETSC_EXTERN PetscErrorCode DMGetVecType(DM,VecType*); 110 PETSC_EXTERN PetscErrorCode DMSetMatType(DM,MatType); 111 PETSC_EXTERN PetscErrorCode DMGetMatType(DM,MatType*); 112 PETSC_EXTERN PetscErrorCode DMSetApplicationContext(DM,void*); 113 PETSC_EXTERN PetscErrorCode DMSetApplicationContextDestroy(DM,PetscErrorCode (*)(void**)); 114 PETSC_EXTERN PetscErrorCode DMGetApplicationContext(DM,void*); 115 PETSC_EXTERN PetscErrorCode DMSetVariableBounds(DM,PetscErrorCode (*)(DM,Vec,Vec)); 116 PETSC_EXTERN PetscErrorCode DMHasVariableBounds(DM,PetscBool *); 117 PETSC_EXTERN PetscErrorCode DMHasColoring(DM,PetscBool *); 118 PETSC_EXTERN PetscErrorCode DMComputeVariableBounds(DM,Vec,Vec); 119 120 PETSC_EXTERN PetscErrorCode DMCreateSubDM(DM, PetscInt, PetscInt[], IS *, DM *); 121 PETSC_EXTERN PetscErrorCode DMCreateFieldDecomposition(DM,PetscInt*,char***,IS**,DM**); 122 PETSC_EXTERN PetscErrorCode DMCreateDomainDecomposition(DM,PetscInt*,char***,IS**,IS**,DM**); 123 PETSC_EXTERN PetscErrorCode DMCreateDomainDecompositionScatters(DM,PetscInt,DM*,VecScatter**,VecScatter**,VecScatter**); 124 125 PETSC_EXTERN PetscErrorCode DMGetRefineLevel(DM,PetscInt*); 126 PETSC_EXTERN PetscErrorCode DMGetCoarsenLevel(DM,PetscInt*); 127 PETSC_EXTERN PetscErrorCode DMFinalizePackage(void); 128 129 PETSC_EXTERN PetscErrorCode VecGetDM(Vec, DM*); 130 PETSC_EXTERN PetscErrorCode VecSetDM(Vec, DM); 131 PETSC_EXTERN PetscErrorCode MatGetDM(Mat, DM*); 132 PETSC_EXTERN PetscErrorCode MatSetDM(Mat, DM); 133 134 typedef struct NLF_DAAD* NLF; 135 136 #define DM_FILE_CLASSID 1211221 137 138 /* FEM support */ 139 PETSC_EXTERN PetscErrorCode DMPrintCellVector(PetscInt, const char [], PetscInt, const PetscScalar []); 140 PETSC_EXTERN PetscErrorCode DMPrintCellMatrix(PetscInt, const char [], PetscInt, PetscInt, const PetscScalar []); 141 PETSC_EXTERN PetscErrorCode DMPrintLocalVec(DM, const char [], PetscReal, Vec); 142 143 PETSC_EXTERN PetscErrorCode DMGetDefaultSection(DM, PetscSection *); 144 PETSC_EXTERN PetscErrorCode DMSetDefaultSection(DM, PetscSection); 145 PETSC_EXTERN PetscErrorCode DMGetDefaultGlobalSection(DM, PetscSection *); 146 PETSC_EXTERN PetscErrorCode DMSetDefaultGlobalSection(DM, PetscSection); 147 PETSC_EXTERN PetscErrorCode DMGetDefaultSF(DM, PetscSF *); 148 PETSC_EXTERN PetscErrorCode DMSetDefaultSF(DM, PetscSF); 149 PETSC_EXTERN PetscErrorCode DMCreateDefaultSF(DM, PetscSection, PetscSection); 150 PETSC_EXTERN PetscErrorCode DMGetPointSF(DM, PetscSF *); 151 PETSC_EXTERN PetscErrorCode DMSetPointSF(DM, PetscSF); 152 153 PETSC_EXTERN PetscErrorCode DMGetOutputDM(DM, DM *); 154 PETSC_EXTERN PetscErrorCode DMGetOutputSequenceNumber(DM, PetscInt *); 155 PETSC_EXTERN PetscErrorCode DMSetOutputSequenceNumber(DM, PetscInt); 156 157 PETSC_EXTERN PetscErrorCode DMGetNumFields(DM, PetscInt *); 158 PETSC_EXTERN PetscErrorCode DMSetNumFields(DM, PetscInt); 159 PETSC_EXTERN PetscErrorCode DMGetField(DM, PetscInt, PetscObject *); 160 PETSC_EXTERN PetscErrorCode DMSetField(DM, PetscInt, PetscObject); 161 162 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; 163 164 struct _DMInterpolationInfo { 165 MPI_Comm comm; 166 PetscInt dim; /*1 The spatial dimension of points */ 167 PetscInt nInput; /* The number of input points */ 168 PetscReal *points; /* The input point coordinates */ 169 PetscInt *cells; /* The cell containing each point */ 170 PetscInt n; /* The number of local points */ 171 Vec coords; /* The point coordinates */ 172 PetscInt dof; /* The number of components to interpolate */ 173 }; 174 typedef struct _DMInterpolationInfo *DMInterpolationInfo; 175 176 PETSC_EXTERN PetscErrorCode DMInterpolationCreate(MPI_Comm, DMInterpolationInfo *); 177 PETSC_EXTERN PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo, PetscInt); 178 PETSC_EXTERN PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo, PetscInt *); 179 PETSC_EXTERN PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo, PetscInt); 180 PETSC_EXTERN PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo, PetscInt *); 181 PETSC_EXTERN PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo, PetscInt, PetscReal[]); 182 PETSC_EXTERN PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo, DM, PetscBool); 183 PETSC_EXTERN PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo, Vec *); 184 PETSC_EXTERN PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo, Vec *); 185 PETSC_EXTERN PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo, Vec *); 186 PETSC_EXTERN PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo, DM, Vec, Vec); 187 PETSC_EXTERN PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *); 188 #endif 189