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 8 PETSC_EXTERN PetscErrorCode DMInitializePackage(const char[]); 9 /*S 10 DM - Abstract PETSc object that manages an abstract grid object and its interactions with the algebraic solvers 11 12 Level: intermediate 13 14 Concepts: grids, grid refinement 15 16 Notes: The DMDACreate() based object and the DMCompositeCreate() based object are examples of DMs 17 18 .seealso: DMCompositeCreate(), DMDACreate(), DMSetType(), DMType 19 S*/ 20 typedef struct _p_DM* DM; 21 22 PETSC_EXTERN PetscClassId DM_CLASSID; 23 24 /*J 25 DMType - String with the name of a PETSc DM or the creation function 26 with an optional dynamic library name, for example 27 http://www.mcs.anl.gov/petsc/lib.a:mydmcreate() 28 29 Level: beginner 30 31 .seealso: DMSetType(), DM 32 J*/ 33 typedef const char* DMType; 34 #define DMDA "da" 35 #define DMADDA "adda" 36 #define DMCOMPOSITE "composite" 37 #define DMSLICED "sliced" 38 #define DMSHELL "shell" 39 #define DMMESH "mesh" 40 #define DMPLEX "complex" 41 #define DMCARTESIAN "cartesian" 42 #define DMREDUNDANT "redundant" 43 #define DMAKKT "akkt" 44 #define DMPATCH "patch" 45 46 PETSC_EXTERN PetscFList DMList; 47 PETSC_EXTERN PetscBool DMRegisterAllCalled; 48 PETSC_EXTERN PetscErrorCode DMCreate(MPI_Comm,DM*); 49 PETSC_EXTERN PetscErrorCode DMSetType(DM, DMType); 50 PETSC_EXTERN PetscErrorCode DMGetType(DM, DMType *); 51 PETSC_EXTERN PetscErrorCode DMRegister(const char[],const char[],const char[],PetscErrorCode (*)(DM)); 52 PETSC_EXTERN PetscErrorCode DMRegisterAll(const char []); 53 PETSC_EXTERN PetscErrorCode DMRegisterDestroy(void); 54 55 56 /*MC 57 DMRegisterDynamic - Adds a new DM component implementation 58 59 Synopsis: 60 PetscErrorCode DMRegisterDynamic(const char *name,const char *path,const char *func_name, PetscErrorCode (*create_func)(DM)) 61 62 Not Collective 63 64 Input Parameters: 65 + name - The name of a new user-defined creation routine 66 . path - The path (either absolute or relative) of the library containing this routine 67 . func_name - The name of routine to create method context 68 - create_func - The creation routine itself 69 70 Notes: 71 DMRegisterDynamic() may be called multiple times to add several user-defined DMs 72 73 If dynamic libraries are used, then the fourth input argument (routine_create) is ignored. 74 75 Sample usage: 76 .vb 77 DMRegisterDynamic("my_da","/home/username/my_lib/lib/libO/solaris/libmy.a", "MyDMCreate", MyDMCreate); 78 .ve 79 80 Then, your DM type can be chosen with the procedural interface via 81 .vb 82 DMCreate(MPI_Comm, DM *); 83 DMSetType(DM,"my_da_name"); 84 .ve 85 or at runtime via the option 86 .vb 87 -da_type my_da_name 88 .ve 89 90 Notes: $PETSC_ARCH occuring in pathname will be replaced with appropriate values. 91 If your function is not being put into a shared library then use DMRegister() instead 92 93 Level: advanced 94 95 .keywords: DM, register 96 .seealso: DMRegisterAll(), DMRegisterDestroy(), DMRegister() 97 M*/ 98 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 99 #define DMRegisterDynamic(a,b,c,d) DMRegister(a,b,c,0) 100 #else 101 #define DMRegisterDynamic(a,b,c,d) DMRegister(a,b,c,d) 102 #endif 103 104 PETSC_EXTERN PetscErrorCode DMView(DM,PetscViewer); 105 PETSC_EXTERN PetscErrorCode DMLoad(DM,PetscViewer); 106 PETSC_EXTERN PetscErrorCode DMDestroy(DM*); 107 PETSC_EXTERN PetscErrorCode DMCreateGlobalVector(DM,Vec*); 108 PETSC_EXTERN PetscErrorCode DMCreateLocalVector(DM,Vec*); 109 PETSC_EXTERN PetscErrorCode DMGetLocalVector(DM,Vec *); 110 PETSC_EXTERN PetscErrorCode DMRestoreLocalVector(DM,Vec *); 111 PETSC_EXTERN PetscErrorCode DMGetGlobalVector(DM,Vec *); 112 PETSC_EXTERN PetscErrorCode DMRestoreGlobalVector(DM,Vec *); 113 PETSC_EXTERN PetscErrorCode DMClearGlobalVectors(DM); 114 PETSC_EXTERN PetscErrorCode DMGetNamedGlobalVector(DM,const char*,Vec*); 115 PETSC_EXTERN PetscErrorCode DMRestoreNamedGlobalVector(DM,const char*,Vec*); 116 PETSC_EXTERN PetscErrorCode DMGetLocalToGlobalMapping(DM,ISLocalToGlobalMapping*); 117 PETSC_EXTERN PetscErrorCode DMGetLocalToGlobalMappingBlock(DM,ISLocalToGlobalMapping*); 118 PETSC_EXTERN PetscErrorCode DMCreateFieldIS(DM,PetscInt*,char***,IS**); 119 PETSC_EXTERN PetscErrorCode DMGetBlockSize(DM,PetscInt*); 120 PETSC_EXTERN PetscErrorCode DMCreateColoring(DM,ISColoringType,MatType,ISColoring*); 121 PETSC_EXTERN PetscErrorCode DMCreateMatrix(DM,MatType,Mat*); 122 PETSC_EXTERN PetscErrorCode DMSetMatrixPreallocateOnly(DM,PetscBool); 123 PETSC_EXTERN PetscErrorCode DMCreateInterpolation(DM,DM,Mat*,Vec*); 124 PETSC_EXTERN PetscErrorCode DMCreateInjection(DM,DM,VecScatter*); 125 PETSC_EXTERN PetscErrorCode DMGetWorkArray(DM,PetscInt,PetscDataType,void*); 126 PETSC_EXTERN PetscErrorCode DMRestoreWorkArray(DM,PetscInt,PetscDataType,void*); 127 PETSC_EXTERN PetscErrorCode DMRefine(DM,MPI_Comm,DM*); 128 PETSC_EXTERN PetscErrorCode DMCoarsen(DM,MPI_Comm,DM*); 129 PETSC_EXTERN PetscErrorCode DMRefineHierarchy(DM,PetscInt,DM[]); 130 PETSC_EXTERN PetscErrorCode DMCoarsenHierarchy(DM,PetscInt,DM[]); 131 PETSC_EXTERN PetscErrorCode DMCoarsenHookAdd(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,Mat,Vec,Mat,DM,void*),void*); 132 PETSC_EXTERN PetscErrorCode DMRefineHookAdd(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,Mat,DM,void*),void*); 133 PETSC_EXTERN PetscErrorCode DMRestrict(DM,Mat,Vec,Mat,DM); 134 PETSC_EXTERN PetscErrorCode DMInterpolate(DM,Mat,DM); 135 PETSC_EXTERN PetscErrorCode DMSetFromOptions(DM); 136 PETSC_EXTERN PetscErrorCode DMViewFromOptions(DM,const char[]); 137 138 PETSC_EXTERN PetscErrorCode DMSetUp(DM); 139 PETSC_EXTERN PetscErrorCode DMCreateInterpolationScale(DM,DM,Mat,Vec*); 140 PETSC_EXTERN PetscErrorCode DMCreateAggregates(DM,DM,Mat*); 141 PETSC_EXTERN PetscErrorCode DMGlobalToLocalHookAdd(DM,PetscErrorCode (*)(DM,Vec,InsertMode,Vec,void*),PetscErrorCode (*)(DM,Vec,InsertMode,Vec,void*),void*); 142 PETSC_EXTERN PetscErrorCode DMGlobalToLocalBegin(DM,Vec,InsertMode,Vec); 143 PETSC_EXTERN PetscErrorCode DMGlobalToLocalEnd(DM,Vec,InsertMode,Vec); 144 PETSC_EXTERN PetscErrorCode DMLocalToGlobalBegin(DM,Vec,InsertMode,Vec); 145 PETSC_EXTERN PetscErrorCode DMLocalToGlobalEnd(DM,Vec,InsertMode,Vec); 146 PETSC_EXTERN PetscErrorCode DMConvert(DM,DMType,DM*); 147 148 PETSC_EXTERN PetscErrorCode DMGetCoordinateDM(DM,DM*); 149 PETSC_EXTERN PetscErrorCode DMGetCoordinates(DM,Vec*); 150 PETSC_EXTERN PetscErrorCode DMSetCoordinates(DM,Vec); 151 PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocal(DM,Vec*); 152 PETSC_EXTERN PetscErrorCode DMSetCoordinatesLocal(DM,Vec); 153 PETSC_EXTERN PetscErrorCode DMLocatePoints(DM,Vec,IS*); 154 155 /* block hook interface */ 156 PETSC_EXTERN PetscErrorCode DMSubDomainHookAdd(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,VecScatter,VecScatter,DM,void*),void*); 157 PETSC_EXTERN PetscErrorCode DMSubDomainRestrict(DM,VecScatter,VecScatter,DM); 158 159 PETSC_EXTERN PetscErrorCode DMSetOptionsPrefix(DM,const char []); 160 PETSC_EXTERN PetscErrorCode DMSetVecType(DM,VecType); 161 PETSC_EXTERN PetscErrorCode DMSetMatType(DM,MatType); 162 PETSC_EXTERN PetscErrorCode DMSetApplicationContext(DM,void*); 163 PETSC_EXTERN PetscErrorCode DMSetApplicationContextDestroy(DM,PetscErrorCode (*)(void**)); 164 PETSC_EXTERN PetscErrorCode DMGetApplicationContext(DM,void*); 165 PETSC_EXTERN PetscErrorCode DMSetVariableBounds(DM,PetscErrorCode (*)(DM,Vec,Vec)); 166 PETSC_EXTERN PetscErrorCode DMHasVariableBounds(DM,PetscBool *); 167 PETSC_EXTERN PetscErrorCode DMHasColoring(DM,PetscBool *); 168 PETSC_EXTERN PetscErrorCode DMComputeVariableBounds(DM,Vec,Vec); 169 170 PETSC_EXTERN PetscErrorCode DMCreateSubDM(DM, PetscInt, PetscInt[], IS *, DM *); 171 PETSC_EXTERN PetscErrorCode DMCreateFieldDecompositionDM(DM,const char*,DM*); 172 PETSC_EXTERN PetscErrorCode DMCreateFieldDecomposition(DM,PetscInt*,char***,IS**,DM**); 173 PETSC_EXTERN PetscErrorCode DMCreateDomainDecompositionDM(DM,const char*,DM*); 174 PETSC_EXTERN PetscErrorCode DMCreateDomainDecomposition(DM,PetscInt*,char***,IS**,IS**,DM**); 175 PETSC_EXTERN PetscErrorCode DMCreateDomainDecompositionScatters(DM,PetscInt,DM*,VecScatter**,VecScatter**,VecScatter**); 176 177 PETSC_EXTERN PetscErrorCode DMGetRefineLevel(DM,PetscInt*); 178 PETSC_EXTERN PetscErrorCode DMGetCoarsenLevel(DM,PetscInt*); 179 PETSC_EXTERN PetscErrorCode DMFinalizePackage(void); 180 181 PETSC_EXTERN PetscErrorCode VecGetDM(Vec, DM*); 182 PETSC_EXTERN PetscErrorCode VecSetDM(Vec, DM); 183 PETSC_EXTERN PetscErrorCode MatGetDM(Mat, DM*); 184 PETSC_EXTERN PetscErrorCode MatSetDM(Mat, DM); 185 186 typedef struct NLF_DAAD* NLF; 187 188 #include <petscbag.h> 189 190 PETSC_EXTERN PetscErrorCode PetscViewerBinaryMatlabOpen(MPI_Comm, const char [], PetscViewer*); 191 PETSC_EXTERN PetscErrorCode PetscViewerBinaryMatlabDestroy(PetscViewer*); 192 PETSC_EXTERN PetscErrorCode PetscViewerBinaryMatlabOutputBag(PetscViewer, const char [], PetscBag); 193 PETSC_EXTERN PetscErrorCode PetscViewerBinaryMatlabOutputVec(PetscViewer, const char [], Vec); 194 PETSC_EXTERN PetscErrorCode PetscViewerBinaryMatlabOutputVecDA(PetscViewer, const char [], Vec, DM); 195 196 #define DM_FILE_CLASSID 1211221 197 198 /* FEM support */ 199 PETSC_EXTERN PetscErrorCode DMPrintCellVector(PetscInt, const char [], PetscInt, const PetscScalar []); 200 PETSC_EXTERN PetscErrorCode DMPrintCellMatrix(PetscInt, const char [], PetscInt, PetscInt, const PetscScalar []); 201 202 PETSC_EXTERN PetscErrorCode DMGetDefaultSection(DM, PetscSection *); 203 PETSC_EXTERN PetscErrorCode DMSetDefaultSection(DM, PetscSection); 204 PETSC_EXTERN PetscErrorCode DMGetDefaultGlobalSection(DM, PetscSection *); 205 PETSC_EXTERN PetscErrorCode DMSetDefaultGlobalSection(DM, PetscSection); 206 PETSC_EXTERN PetscErrorCode DMGetDefaultSF(DM, PetscSF *); 207 PETSC_EXTERN PetscErrorCode DMSetDefaultSF(DM, PetscSF); 208 PETSC_EXTERN PetscErrorCode DMCreateDefaultSF(DM, PetscSection, PetscSection); 209 PETSC_EXTERN PetscErrorCode DMGetPointSF(DM, PetscSF *); 210 PETSC_EXTERN PetscErrorCode DMSetPointSF(DM, PetscSF); 211 212 PETSC_EXTERN PetscErrorCode DMGetNumFields(DM, PetscInt *); 213 PETSC_EXTERN PetscErrorCode DMSetNumFields(DM, PetscInt); 214 PETSC_EXTERN PetscErrorCode DMGetField(DM, PetscInt, PetscObject *); 215 216 typedef struct { 217 PetscInt numQuadPoints; /* The number of quadrature points on an element */ 218 const PetscReal *quadPoints; /* The quadrature point coordinates */ 219 const PetscReal *quadWeights; /* The quadrature weights */ 220 PetscInt numBasisFuncs; /* The number of finite element basis functions on an element */ 221 PetscInt numComponents; /* The number of components for each basis function */ 222 const PetscReal *basis; /* The basis functions tabulated at the quadrature points */ 223 const PetscReal *basisDer; /* The basis function derivatives tabulated at the quadrature points */ 224 } PetscQuadrature; 225 226 typedef struct { 227 PetscQuadrature *quad; 228 void (**f0Funcs)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The f_0 functions for each field */ 229 void (**f1Funcs)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The f_1 functions for each field */ 230 void (**g0Funcs)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The g_0 functions for each field pair */ 231 void (**g1Funcs)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The g_1 functions for each field pair */ 232 void (**g2Funcs)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The g_2 functions for each field pair */ 233 void (**g3Funcs)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The g_3 functions for each field pair */ 234 PetscScalar (**bcFuncs)(const PetscReal x[]); /* The boundary condition function for each field component */ 235 } PetscFEM; 236 237 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; 238 239 struct _DMInterpolationInfo { 240 MPI_Comm comm; 241 PetscInt dim; /*1 The spatial dimension of points */ 242 PetscInt nInput; /* The number of input points */ 243 PetscReal *points; /* The input point coordinates */ 244 PetscInt *cells; /* The cell containing each point */ 245 PetscInt n; /* The number of local points */ 246 Vec coords; /* The point coordinates */ 247 PetscInt dof; /* The number of components to interpolate */ 248 }; 249 typedef struct _DMInterpolationInfo *DMInterpolationInfo; 250 251 PetscErrorCode DMInterpolationCreate(MPI_Comm, DMInterpolationInfo *); 252 PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo, PetscInt); 253 PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo, PetscInt *); 254 PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo, PetscInt); 255 PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo, PetscInt *); 256 PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo, PetscInt, PetscReal[]); 257 PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo, DM, PetscBool); 258 PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo, Vec *); 259 PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo, Vec *); 260 PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo, Vec *); 261 PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo, DM, Vec, Vec); 262 PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *); 263 #endif 264