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