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