1 /* 2 Objects to manage the interactions between the mesh data structures and the algebraic objects 3 */ 4 #if !defined(__PETSCDA_H) 5 #define __PETSCDA_H 6 #include "petscmat.h" 7 #include "petscao.h" 8 PETSC_EXTERN_CXX_BEGIN 9 10 extern PetscErrorCode DMInitializePackage(const char[]); 11 /*S 12 DM - Abstract PETSc object that manages an abstract grid object 13 14 Level: intermediate 15 16 Concepts: grids, grid refinement 17 18 Notes: The DMDACreate() based object and the DMCompositeCreate() based object are examples of DMs 19 20 Though the DM objects require the petscsnes.h include files the DM library is 21 NOT dependent on the SNES or KSP library. In fact, the KSP and SNES libraries depend on 22 DM. (This is not great design, but not trivial to fix). 23 24 .seealso: DMCompositeCreate(), DMDACreate() 25 S*/ 26 typedef struct _p_DM* DM; 27 28 /*E 29 DMDAStencilType - Determines if the stencil extends only along the coordinate directions, or also 30 to the northeast, northwest etc 31 32 Level: beginner 33 34 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDACreate() 35 E*/ 36 typedef enum { DMDA_STENCIL_STAR,DMDA_STENCIL_BOX } DMDAStencilType; 37 38 /*MC 39 DMDA_STENCIL_STAR - "Star"-type stencil. In logical grid coordinates, only (i,j,k), (i+s,j,k), (i,j+s,k), 40 (i,j,k+s) are in the stencil NOT, for example, (i+s,j+s,k) 41 42 Level: beginner 43 44 .seealso: DMDA_STENCIL_BOX, DMDAStencilType 45 M*/ 46 47 /*MC 48 DMDA_STENCIL_BOX - "Box"-type stencil. In logical grid coordinates, any of (i,j,k), (i+s,j+r,k+t) may 49 be in the stencil. 50 51 Level: beginner 52 53 .seealso: DMDA_STENCIL_STAR, DMDAStencilType 54 M*/ 55 56 /*E 57 DMDABoundaryType - Describes the choice for fill of ghost cells on domain boundaries. 58 59 Level: beginner 60 61 A boundary may be of type DMDA_BOUNDARY_NONE (no ghost nodes), DMDA_BOUNDARY_GHOST (ghost nodes 62 exist but aren't filled), DMDA_BOUNDARY_MIRROR (not yet implemented), or DMDA_BOUNDARY_PERIODIC 63 (ghost nodes filled by the opposite edge of the domain). 64 65 .seealso: DMDASetBoundaryType(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDACreate() 66 E*/ 67 typedef enum { DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_GHOSTED, DMDA_BOUNDARY_MIRROR, DMDA_BOUNDARY_PERIODIC } DMDABoundaryType; 68 69 extern const char *DMDABoundaryTypes[]; 70 71 /*E 72 DMDAInterpolationType - Defines the type of interpolation that will be returned by 73 DMGetInterpolation. 74 75 Level: beginner 76 77 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMGetInterpolation(), DMDASetInterpolationType(), DMDACreate() 78 E*/ 79 typedef enum { DMDA_Q0, DMDA_Q1 } DMDAInterpolationType; 80 81 extern PetscErrorCode DMDASetInterpolationType(DM,DMDAInterpolationType); 82 83 /*E 84 DMDAElementType - Defines the type of elements that will be returned by 85 DMGetElements() 86 87 Level: beginner 88 89 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMGetInterpolation(), DMDASetInterpolationType(), 90 DMDASetElementType(), DMGetElements(), DMRestoreElements(), DMDACreate() 91 E*/ 92 typedef enum { DMDA_ELEMENT_P1, DMDA_ELEMENT_Q1 } DMDAElementType; 93 94 extern PetscErrorCode DMDASetElementType(DM,DMDAElementType); 95 extern PetscErrorCode DMDAGetElementType(DM,DMDAElementType*); 96 97 typedef enum { DMDA_X,DMDA_Y,DMDA_Z } DMDADirection; 98 99 extern PetscClassId DM_CLASSID; 100 101 #define MATSEQUSFFT "sequsfft" 102 103 extern PetscErrorCode DMDACreate(MPI_Comm,DM*); 104 extern PetscErrorCode DMDASetDim(DM,PetscInt); 105 extern PetscErrorCode DMDASetSizes(DM,PetscInt,PetscInt,PetscInt); 106 extern PetscErrorCode DMDACreate1d(MPI_Comm,DMDABoundaryType,PetscInt,PetscInt,PetscInt,const PetscInt[],DM *); 107 extern PetscErrorCode DMDACreate2d(MPI_Comm,DMDABoundaryType,DMDABoundaryType,DMDAStencilType,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],DM*); 108 extern PetscErrorCode DMDACreate3d(MPI_Comm,DMDABoundaryType,DMDABoundaryType,DMDABoundaryType,DMDAStencilType,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscInt[],DM*); 109 extern PetscErrorCode DMSetOptionsPrefix(DM,const char []); 110 extern PetscErrorCode DMSetVecType(DM,const VecType); 111 112 extern PetscErrorCode DMDAGlobalToNaturalBegin(DM,Vec,InsertMode,Vec); 113 extern PetscErrorCode DMDAGlobalToNaturalEnd(DM,Vec,InsertMode,Vec); 114 extern PetscErrorCode DMDANaturalToGlobalBegin(DM,Vec,InsertMode,Vec); 115 extern PetscErrorCode DMDANaturalToGlobalEnd(DM,Vec,InsertMode,Vec); 116 extern PetscErrorCode DMDALocalToLocalBegin(DM,Vec,InsertMode,Vec); 117 extern PetscErrorCode DMDALocalToLocalEnd(DM,Vec,InsertMode,Vec); 118 extern PetscErrorCode DMDACreateNaturalVector(DM,Vec *); 119 120 extern PetscErrorCode DMDALoad(PetscViewer,PetscInt,PetscInt,PetscInt,DM *); 121 extern PetscErrorCode DMDAGetCorners(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*); 122 extern PetscErrorCode DMDAGetGhostCorners(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*); 123 extern PetscErrorCode DMDAGetInfo(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,DMDABoundaryType*,DMDABoundaryType*,DMDABoundaryType*,DMDAStencilType*); 124 extern PetscErrorCode DMDAGetProcessorSubset(DM,DMDADirection,PetscInt,MPI_Comm*); 125 extern PetscErrorCode DMDAGetProcessorSubsets(DM,DMDADirection,MPI_Comm*); 126 127 extern PetscErrorCode DMDAGlobalToNaturalAllCreate(DM,VecScatter*); 128 extern PetscErrorCode DMDANaturalAllToGlobalCreate(DM,VecScatter*); 129 130 extern PetscErrorCode DMDAGetGlobalIndices(DM,PetscInt*,PetscInt**); 131 132 extern PetscErrorCode DMDAGetScatter(DM,VecScatter*,VecScatter*,VecScatter*); 133 extern PetscErrorCode DMDAGetNeighbors(DM,const PetscMPIInt**); 134 135 extern PetscErrorCode DMDAGetAO(DM,AO*); 136 extern PetscErrorCode DMDASetCoordinates(DM,Vec); 137 extern PetscErrorCode DMDASetGhostedCoordinates(DM,Vec); 138 extern PetscErrorCode DMDAGetCoordinates(DM,Vec *); 139 extern PetscErrorCode DMDAGetGhostedCoordinates(DM,Vec *); 140 extern PetscErrorCode DMDAGetCoordinateDA(DM,DM *); 141 extern PetscErrorCode DMDASetUniformCoordinates(DM,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal); 142 extern PetscErrorCode DMDAGetBoundingBox(DM,PetscReal[],PetscReal[]); 143 extern PetscErrorCode DMDAGetLocalBoundingBox(DM,PetscReal[],PetscReal[]); 144 extern PetscErrorCode DMDASetFieldName(DM,PetscInt,const char[]); 145 extern PetscErrorCode DMDAGetFieldName(DM,PetscInt,const char**); 146 147 extern PetscErrorCode DMDASetBoundaryType(DM,DMDABoundaryType,DMDABoundaryType,DMDABoundaryType); 148 extern PetscErrorCode DMDASetDof(DM, int); 149 extern PetscErrorCode DMDASetStencilWidth(DM, PetscInt); 150 extern PetscErrorCode DMDASetOwnershipRanges(DM,const PetscInt[],const PetscInt[],const PetscInt[]); 151 extern PetscErrorCode DMDAGetOwnershipRanges(DM,const PetscInt**,const PetscInt**,const PetscInt**); 152 extern PetscErrorCode DMDASetNumProcs(DM, PetscInt, PetscInt, PetscInt); 153 extern PetscErrorCode DMDASetStencilType(DM, DMDAStencilType); 154 155 extern PetscErrorCode DMDAVecGetArray(DM,Vec,void *); 156 extern PetscErrorCode DMDAVecRestoreArray(DM,Vec,void *); 157 158 extern PetscErrorCode DMDAVecGetArrayDOF(DM,Vec,void *); 159 extern PetscErrorCode DMDAVecRestoreArrayDOF(DM,Vec,void *); 160 161 extern PetscErrorCode DMDASplitComm2d(MPI_Comm,PetscInt,PetscInt,PetscInt,MPI_Comm*); 162 163 /*E 164 DMType - String with the name of a PETSc DM or the creation function 165 with an optional dynamic library name, for example 166 http://www.mcs.anl.gov/petsc/lib.a:myveccreate() 167 168 Level: beginner 169 170 .seealso: DMSetType(), DM 171 E*/ 172 173 #define DMType char* 174 #define DMDA "da" 175 #define DMADDA "adda" 176 #define DMCOMPOSITE "composite" 177 #define DMSLICED "sliced" 178 #define DMMESH "mesh" 179 #define DMCARTESIAN "cartesian" 180 181 extern PetscFList DMList; 182 extern PetscBool DMRegisterAllCalled; 183 extern PetscErrorCode DMCreate(MPI_Comm,DM*); 184 extern PetscErrorCode DMSetType(DM, const DMType); 185 extern PetscErrorCode DMGetType(DM, const DMType *); 186 extern PetscErrorCode DMRegister(const char[],const char[],const char[],PetscErrorCode (*)(DM)); 187 extern PetscErrorCode DMRegisterAll(const char []); 188 extern PetscErrorCode DMRegisterDestroy(void); 189 190 /*MC 191 DMRegisterDynamic - Adds a new DM component implementation 192 193 Synopsis: 194 PetscErrorCode DMRegisterDynamic(const char *name,const char *path,const char *func_name, PetscErrorCode (*create_func)(DM)) 195 196 Not Collective 197 198 Input Parameters: 199 + name - The name of a new user-defined creation routine 200 . path - The path (either absolute or relative) of the library containing this routine 201 . func_name - The name of routine to create method context 202 - create_func - The creation routine itself 203 204 Notes: 205 DMRegisterDynamic() may be called multiple times to add several user-defined DMs 206 207 If dynamic libraries are used, then the fourth input argument (routine_create) is ignored. 208 209 Sample usage: 210 .vb 211 DMRegisterDynamic("my_da","/home/username/my_lib/lib/libO/solaris/libmy.a", "MyDMCreate", MyDMCreate); 212 .ve 213 214 Then, your DM type can be chosen with the procedural interface via 215 .vb 216 DMCreate(MPI_Comm, DM *); 217 DMSetType(DM,"my_da_name"); 218 .ve 219 or at runtime via the option 220 .vb 221 -da_type my_da_name 222 .ve 223 224 Notes: $PETSC_ARCH occuring in pathname will be replaced with appropriate values. 225 If your function is not being put into a shared library then use DMRegister() instead 226 227 Level: advanced 228 229 .keywords: DM, register 230 .seealso: DMRegisterAll(), DMRegisterDestroy(), DMRegister() 231 M*/ 232 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 233 #define DMRegisterDynamic(a,b,c,d) DMRegister(a,b,c,0) 234 #else 235 #define DMRegisterDynamic(a,b,c,d) DMRegister(a,b,c,d) 236 #endif 237 238 extern PetscErrorCode MatRegisterDAAD(void); 239 extern PetscErrorCode MatCreateDAAD(DM,Mat*); 240 extern PetscErrorCode MatCreateSeqUSFFT(Vec, DM,Mat*); 241 242 /*S 243 DMDALocalInfo - C struct that contains information about a structured grid and a processors logical 244 location in it. 245 246 Level: beginner 247 248 Concepts: distributed array 249 250 Developer note: Then entries in this struct are int instead of PetscInt so that the elements may 251 be extracted in Fortran as if from an integer array 252 253 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DM, DMDAGetLocalInfo(), DMDAGetInfo() 254 S*/ 255 typedef struct { 256 PetscInt dim,dof,sw; 257 PetscInt mx,my,mz; /* global number of grid points in each direction */ 258 PetscInt xs,ys,zs; /* starting pointd of this processor, excluding ghosts */ 259 PetscInt xm,ym,zm; /* number of grid points on this processor, excluding ghosts */ 260 PetscInt gxs,gys,gzs; /* starting point of this processor including ghosts */ 261 PetscInt gxm,gym,gzm; /* number of grid points on this processor including ghosts */ 262 DMDABoundaryType bx,by,bz; /* type of ghost nodes at boundary */ 263 DMDAStencilType st; 264 DM da; 265 } DMDALocalInfo; 266 267 /*MC 268 DMDAForEachPointBegin2d - Starts a loop over the local part of a two dimensional DMDA 269 270 Synopsis: 271 void DMDAForEachPointBegin2d(DALocalInfo *info,PetscInt i,PetscInt j); 272 273 Not Collective 274 275 Level: intermediate 276 277 .seealso: DMDAForEachPointEnd2d(), DMDAVecGetArray() 278 M*/ 279 #define DMDAForEachPointBegin2d(info,i,j) {\ 280 PetscInt _xints = info->xs,_xinte = info->xs+info->xm,_yints = info->ys,_yinte = info->ys+info->ym;\ 281 for (j=_yints; j<_yinte; j++) {\ 282 for (i=_xints; i<_xinte; i++) {\ 283 284 /*MC 285 DMDAForEachPointEnd2d - Ends a loop over the local part of a two dimensional DMDA 286 287 Synopsis: 288 void DMDAForEachPointEnd2d; 289 290 Not Collective 291 292 Level: intermediate 293 294 .seealso: DMDAForEachPointBegin2d(), DMDAVecGetArray() 295 M*/ 296 #define DMDAForEachPointEnd2d }}} 297 298 /*MC 299 DMDACoor2d - Structure for holding 2d (x and y) coordinates. 300 301 Level: intermediate 302 303 Sample Usage: 304 DMDACoor2d **coors; 305 Vec vcoors; 306 DM cda; 307 308 DMDAGetCoordinates(da,&vcoors); 309 DMDAGetCoordinateDA(da,&cda); 310 DMDAVecGetArray(cda,vcoors,&coors); 311 DMDAGetCorners(cda,&mstart,&nstart,0,&m,&n,0) 312 for (i=mstart; i<mstart+m; i++) { 313 for (j=nstart; j<nstart+n; j++) { 314 x = coors[j][i].x; 315 y = coors[j][i].y; 316 ...... 317 } 318 } 319 DMDAVecRestoreArray(dac,vcoors,&coors); 320 321 .seealso: DMDACoor3d, DMDAForEachPointBegin(), DMDAGetCoordinateDA(), DMDAGetCoordinates(), DMDAGetGhostCoordinates() 322 M*/ 323 typedef struct {PetscScalar x,y;} DMDACoor2d; 324 325 /*MC 326 DMDACoor3d - Structure for holding 3d (x, y and z) coordinates. 327 328 Level: intermediate 329 330 Sample Usage: 331 DMDACoor3d ***coors; 332 Vec vcoors; 333 DM cda; 334 335 DMDAGetCoordinates(da,&vcoors); 336 DMDAGetCoordinateDA(da,&cda); 337 DMDAVecGetArray(cda,vcoors,&coors); 338 DMDAGetCorners(cda,&mstart,&nstart,&pstart,&m,&n,&p) 339 for (i=mstart; i<mstart+m; i++) { 340 for (j=nstart; j<nstart+n; j++) { 341 for (k=pstart; k<pstart+p; k++) { 342 x = coors[k][j][i].x; 343 y = coors[k][j][i].y; 344 z = coors[k][j][i].z; 345 ...... 346 } 347 } 348 DMDAVecRestoreArray(dac,vcoors,&coors); 349 350 .seealso: DMDACoor2d, DMDAForEachPointBegin(), DMDAGetCoordinateDA(), DMDAGetCoordinates(), DMDAGetGhostCoordinates() 351 M*/ 352 typedef struct {PetscScalar x,y,z;} DMDACoor3d; 353 354 extern PetscErrorCode DMDAGetLocalInfo(DM,DMDALocalInfo*); 355 typedef PetscErrorCode (*DMDALocalFunction1)(DMDALocalInfo*,void*,void*,void*); 356 extern PetscErrorCode DMDAFormFunctionLocal(DM, DMDALocalFunction1, Vec, Vec, void *); 357 extern PetscErrorCode DMDAFormFunctionLocalGhost(DM, DMDALocalFunction1, Vec, Vec, void *); 358 extern PetscErrorCode DMDAFormJacobianLocal(DM, DMDALocalFunction1, Vec, Mat, void *); 359 extern PetscErrorCode DMDAFormFunction1(DM,Vec,Vec,void*); 360 extern PetscErrorCode DMDAFormFunction(DM,PetscErrorCode (*)(void),Vec,Vec,void*); 361 extern PetscErrorCode DMDAFormFunctioni1(DM,PetscInt,Vec,PetscScalar*,void*); 362 extern PetscErrorCode DMDAFormFunctionib1(DM,PetscInt,Vec,PetscScalar*,void*); 363 extern PetscErrorCode DMDAComputeJacobian1WithAdic(DM,Vec,Mat,void*); 364 extern PetscErrorCode DMDAComputeJacobian1WithAdifor(DM,Vec,Mat,void*); 365 extern PetscErrorCode DMDAMultiplyByJacobian1WithAdic(DM,Vec,Vec,Vec,void*); 366 extern PetscErrorCode DMDAMultiplyByJacobian1WithAdifor(DM,Vec,Vec,Vec,void*); 367 extern PetscErrorCode DMDAMultiplyByJacobian1WithAD(DM,Vec,Vec,Vec,void*); 368 extern PetscErrorCode DMDAComputeJacobian1(DM,Vec,Mat,void*); 369 extern PetscErrorCode DMDAGetLocalFunction(DM,DMDALocalFunction1*); 370 extern PetscErrorCode DMDASetLocalFunction(DM,DMDALocalFunction1); 371 extern PetscErrorCode DMDASetLocalFunctioni(DM,PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,PetscScalar*,void*)); 372 extern PetscErrorCode DMDASetLocalFunctionib(DM,PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,PetscScalar*,void*)); 373 extern PetscErrorCode DMDAGetLocalJacobian(DM,DMDALocalFunction1*); 374 extern PetscErrorCode DMDASetLocalJacobian(DM,DMDALocalFunction1); 375 extern PetscErrorCode DMDASetLocalAdicFunction_Private(DM,DMDALocalFunction1); 376 377 extern PetscErrorCode MatSetDA(Mat,DM); 378 379 /*MC 380 DMDASetLocalAdicFunction - Caches in a DM a local function computed by ADIC/ADIFOR 381 382 Synopsis: 383 PetscErrorCode DMDASetLocalAdicFunction(DM da,DMDALocalFunction1 ad_lf) 384 385 Logically Collective on DM 386 387 Input Parameter: 388 + da - initial distributed array 389 - ad_lf - the local function as computed by ADIC/ADIFOR 390 391 Level: intermediate 392 393 .keywords: distributed array, refine 394 395 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDAGetLocalFunction(), DMDASetLocalFunction(), 396 DMDASetLocalJacobian() 397 M*/ 398 #if defined(PETSC_HAVE_ADIC) 399 # define DMDASetLocalAdicFunction(a,d) DMDASetLocalAdicFunction_Private(a,(DMDALocalFunction1)d) 400 #else 401 # define DMDASetLocalAdicFunction(a,d) DMDASetLocalAdicFunction_Private(a,0) 402 #endif 403 404 extern PetscErrorCode DMDASetLocalAdicMFFunction_Private(DM,DMDALocalFunction1); 405 #if defined(PETSC_HAVE_ADIC) 406 # define DMDASetLocalAdicMFFunction(a,d) DMDASetLocalAdicMFFunction_Private(a,(DMDALocalFunction1)d) 407 #else 408 # define DMDASetLocalAdicMFFunction(a,d) DMDASetLocalAdicMFFunction_Private(a,0) 409 #endif 410 extern PetscErrorCode DMDASetLocalAdicFunctioni_Private(DM,PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*)); 411 #if defined(PETSC_HAVE_ADIC) 412 # define DMDASetLocalAdicFunctioni(a,d) DMDASetLocalAdicFunctioni_Private(a,(PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*))d) 413 #else 414 # define DMDASetLocalAdicFunctioni(a,d) DMDASetLocalAdicFunctioni_Private(a,0) 415 #endif 416 extern PetscErrorCode DMDASetLocalAdicMFFunctioni_Private(DM,PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*)); 417 #if defined(PETSC_HAVE_ADIC) 418 # define DMDASetLocalAdicMFFunctioni(a,d) DMDASetLocalAdicMFFunctioni_Private(a,(PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*))d) 419 #else 420 # define DMDASetLocalAdicMFFunctioni(a,d) DMDASetLocalAdicMFFunctioni_Private(a,0) 421 #endif 422 423 extern PetscErrorCode DMDASetLocalAdicFunctionib_Private(DM,PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*)); 424 #if defined(PETSC_HAVE_ADIC) 425 # define DMDASetLocalAdicFunctionib(a,d) DMDASetLocalAdicFunctionib_Private(a,(PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*))d) 426 #else 427 # define DMDASetLocalAdicFunctionib(a,d) DMDASetLocalAdicFunctionib_Private(a,0) 428 #endif 429 extern PetscErrorCode DMDASetLocalAdicMFFunctionib_Private(DM,PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*)); 430 #if defined(PETSC_HAVE_ADIC) 431 # define DMDASetLocalAdicMFFunctionib(a,d) DMDASetLocalAdicMFFunctionib_Private(a,(PetscErrorCode (*)(DMDALocalInfo*,MatStencil*,void*,void*,void*))d) 432 #else 433 # define DMDASetLocalAdicMFFunctionib(a,d) DMDASetLocalAdicMFFunctionib_Private(a,0) 434 #endif 435 436 extern PetscErrorCode DMDAFormFunctioniTest1(DM,void*); 437 438 extern PetscErrorCode DMView(DM,PetscViewer); 439 extern PetscErrorCode DMDestroy(DM); 440 extern PetscErrorCode DMCreateGlobalVector(DM,Vec*); 441 extern PetscErrorCode DMCreateLocalVector(DM,Vec*); 442 extern PetscErrorCode DMGetLocalVector(DM,Vec *); 443 extern PetscErrorCode DMRestoreLocalVector(DM,Vec *); 444 extern PetscErrorCode DMGetGlobalVector(DM,Vec *); 445 extern PetscErrorCode DMRestoreGlobalVector(DM,Vec *); 446 extern PetscErrorCode DMGetLocalToGlobalMapping(DM,ISLocalToGlobalMapping*); 447 extern PetscErrorCode DMGetLocalToGlobalMappingBlock(DM,ISLocalToGlobalMapping*); 448 extern PetscErrorCode DMGetBlockSize(DM,PetscInt*); 449 extern PetscErrorCode DMGetColoring(DM,ISColoringType,const MatType,ISColoring*); 450 extern PetscErrorCode DMGetMatrix(DM, const MatType,Mat*); 451 extern PetscErrorCode DMGetInterpolation(DM,DM,Mat*,Vec*); 452 extern PetscErrorCode DMGetInjection(DM,DM,VecScatter*); 453 extern PetscErrorCode DMRefine(DM,MPI_Comm,DM*); 454 extern PetscErrorCode DMCoarsen(DM,MPI_Comm,DM*); 455 extern PetscErrorCode DMRefineHierarchy(DM,PetscInt,DM[]); 456 extern PetscErrorCode DMCoarsenHierarchy(DM,PetscInt,DM[]); 457 extern PetscErrorCode DMSetFromOptions(DM); 458 extern PetscErrorCode DMSetUp(DM); 459 extern PetscErrorCode DMGetInterpolationScale(DM,DM,Mat,Vec*); 460 extern PetscErrorCode DMGetAggregates(DM,DM,Mat*); 461 extern PetscErrorCode DMGlobalToLocalBegin(DM,Vec,InsertMode,Vec); 462 extern PetscErrorCode DMGlobalToLocalEnd(DM,Vec,InsertMode,Vec); 463 extern PetscErrorCode DMLocalToGlobalBegin(DM,Vec,InsertMode,Vec); 464 extern PetscErrorCode DMLocalToGlobalEnd(DM,Vec,InsertMode,Vec); 465 extern PetscErrorCode DMGetElements(DM,PetscInt *,PetscInt *,const PetscInt*[]); 466 extern PetscErrorCode DMRestoreElements(DM,PetscInt *,PetscInt *,const PetscInt*[]); 467 468 extern PetscErrorCode DMSetContext(DM,void*); 469 extern PetscErrorCode DMGetContext(DM,void**); 470 extern PetscErrorCode DMSetInitialGuess(DM,PetscErrorCode (*)(DM,Vec)); 471 extern PetscErrorCode DMSetFunction(DM,PetscErrorCode (*)(DM,Vec,Vec)); 472 extern PetscErrorCode DMSetJacobian(DM,PetscErrorCode (*)(DM,Vec,Mat,Mat,MatStructure *)); 473 extern PetscErrorCode DMHasInitialGuess(DM,PetscBool *); 474 extern PetscErrorCode DMHasFunction(DM,PetscBool *); 475 extern PetscErrorCode DMHasJacobian(DM,PetscBool *); 476 extern PetscErrorCode DMComputeInitialGuess(DM,Vec); 477 extern PetscErrorCode DMComputeFunction(DM,Vec,Vec); 478 extern PetscErrorCode DMComputeJacobian(DM,Vec,Mat,Mat,MatStructure *); 479 extern PetscErrorCode DMComputeJacobianDefault(DM,Vec,Mat,Mat,MatStructure *); 480 extern PetscErrorCode DMFinalizePackage(void); 481 482 extern PetscErrorCode DMDASetGetMatrix(DM,PetscErrorCode (*)(DM, const MatType,Mat *)); 483 extern PetscErrorCode DMDASetBlockFills(DM,PetscInt*,PetscInt*); 484 extern PetscErrorCode DMDASetMatPreallocateOnly(DM,PetscBool ); 485 extern PetscErrorCode DMDASetRefinementFactor(DM,PetscInt,PetscInt,PetscInt); 486 extern PetscErrorCode DMDAGetRefinementFactor(DM,PetscInt*,PetscInt*,PetscInt*); 487 488 extern PetscErrorCode DMDAGetAdicArray(DM,PetscBool ,void*,void*,PetscInt*); 489 extern PetscErrorCode DMDARestoreAdicArray(DM,PetscBool ,void*,void*,PetscInt*); 490 extern PetscErrorCode DMDAGetAdicMFArray(DM,PetscBool ,void*,void*,PetscInt*); 491 extern PetscErrorCode DMDAGetAdicMFArray4(DM,PetscBool ,void*,void*,PetscInt*); 492 extern PetscErrorCode DMDAGetAdicMFArray9(DM,PetscBool ,void*,void*,PetscInt*); 493 extern PetscErrorCode DMDAGetAdicMFArrayb(DM,PetscBool ,void*,void*,PetscInt*); 494 extern PetscErrorCode DMDARestoreAdicMFArray(DM,PetscBool ,void*,void*,PetscInt*); 495 extern PetscErrorCode DMDAGetArray(DM,PetscBool ,void*); 496 extern PetscErrorCode DMDARestoreArray(DM,PetscBool ,void*); 497 extern PetscErrorCode ad_DAGetArray(DM,PetscBool ,void*); 498 extern PetscErrorCode ad_DARestoreArray(DM,PetscBool ,void*); 499 extern PetscErrorCode admf_DAGetArray(DM,PetscBool ,void*); 500 extern PetscErrorCode admf_DARestoreArray(DM,PetscBool ,void*); 501 502 #include "petscpf.h" 503 extern PetscErrorCode DMDACreatePF(DM,PF*); 504 505 extern PetscErrorCode DMCompositeCreate(MPI_Comm,DM*); 506 extern PetscErrorCode DMCompositeAddArray(DM,PetscMPIInt,PetscInt); 507 extern PetscErrorCode DMCompositeAddDM(DM,DM); 508 extern PetscErrorCode DMCompositeSetCoupling(DM,PetscErrorCode (*)(DM,Mat,PetscInt*,PetscInt*,PetscInt,PetscInt,PetscInt,PetscInt)); 509 extern PetscErrorCode DMCompositeSetContext(DM,void*); 510 extern PetscErrorCode DMCompositeGetContext(DM,void**); 511 extern PetscErrorCode DMCompositeAddVecScatter(DM,VecScatter); 512 extern PetscErrorCode DMCompositeScatter(DM,Vec,...); 513 extern PetscErrorCode DMCompositeGather(DM,Vec,InsertMode,...); 514 extern PetscErrorCode DMCompositeGetAccess(DM,Vec,...); 515 extern PetscErrorCode DMCompositeGetNumberDM(DM,PetscInt*); 516 extern PetscErrorCode DMCompositeRestoreAccess(DM,Vec,...); 517 extern PetscErrorCode DMCompositeGetLocalVectors(DM,...); 518 extern PetscErrorCode DMCompositeGetEntries(DM,...); 519 extern PetscErrorCode DMCompositeRestoreLocalVectors(DM,...); 520 extern PetscErrorCode DMCompositeGetGlobalISs(DM,IS*[]); 521 extern PetscErrorCode DMCompositeGetLocalISs(DM,IS**); 522 extern PetscErrorCode DMCompositeGetISLocalToGlobalMappings(DM,ISLocalToGlobalMapping**); 523 524 extern PetscErrorCode DMSlicedCreate(MPI_Comm,DM*); 525 extern PetscErrorCode DMSlicedGetGlobalIndices(DM,PetscInt*[]); 526 extern PetscErrorCode DMSlicedSetPreallocation(DM,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 527 extern PetscErrorCode DMSlicedSetBlockFills(DM,const PetscInt*,const PetscInt*); 528 extern PetscErrorCode DMSlicedSetGhosts(DM,PetscInt,PetscInt,PetscInt,const PetscInt[]); 529 530 531 typedef struct NLF_DAAD* NLF; 532 533 #include <petscbag.h> 534 535 extern PetscErrorCode PetscViewerBinaryMatlabOpen(MPI_Comm, const char [], PetscViewer*); 536 extern PetscErrorCode PetscViewerBinaryMatlabDestroy(PetscViewer); 537 extern PetscErrorCode PetscViewerBinaryMatlabOutputBag(PetscViewer, const char [], PetscBag); 538 extern PetscErrorCode PetscViewerBinaryMatlabOutputVec(PetscViewer, const char [], Vec); 539 extern PetscErrorCode PetscViewerBinaryMatlabOutputVecDA(PetscViewer, const char [], Vec, DM); 540 541 542 PetscErrorCode DMADDACreate(MPI_Comm,PetscInt,PetscInt*,PetscInt*,PetscInt,PetscBool *,DM*); 543 PetscErrorCode DMADDASetParameters(DM,PetscInt,PetscInt*,PetscInt*,PetscInt,PetscBool*); 544 PetscErrorCode DMADDASetRefinement(DM, PetscInt *,PetscInt); 545 PetscErrorCode DMADDAGetCorners(DM, PetscInt **, PetscInt **); 546 PetscErrorCode DMADDAGetGhostCorners(DM, PetscInt **, PetscInt **); 547 PetscErrorCode DMADDAGetMatrixNS(DM, DM, const MatType , Mat *); 548 549 /* functions to set values in vectors and matrices */ 550 struct _ADDAIdx_s { 551 PetscInt *x; /* the coordinates, user has to make sure it is the correct size! */ 552 PetscInt d; /* indexes the dof */ 553 }; 554 typedef struct _ADDAIdx_s ADDAIdx; 555 556 PetscErrorCode DMADDAMatSetValues(Mat, DM, PetscInt, const ADDAIdx[], DM, PetscInt, const ADDAIdx[], const PetscScalar[], InsertMode); 557 PetscBool ADDAHCiterStartup(const PetscInt, const PetscInt *const, const PetscInt *const, PetscInt *const); 558 PetscBool ADDAHCiter(const PetscInt, const PetscInt *const, const PetscInt *const, PetscInt *const); 559 560 /* Mesh Section */ 561 #if defined(PETSC_HAVE_SIEVE) && defined(__cplusplus) 562 563 #include <sieve/Mesh.hh> 564 #include <sieve/CartesianSieve.hh> 565 #include <sieve/Distribution.hh> 566 #include <sieve/Generator.hh> 567 568 extern PetscErrorCode DMMeshCreate(MPI_Comm, DM*); 569 extern PetscErrorCode DMMeshGetMesh(DM, ALE::Obj<PETSC_MESH_TYPE>&); 570 extern PetscErrorCode DMMeshSetMesh(DM, const ALE::Obj<PETSC_MESH_TYPE>&); 571 extern PetscErrorCode DMMeshGetGlobalScatter(DM, VecScatter *); 572 extern PetscErrorCode DMMeshFinalize(); 573 574 extern PetscErrorCode DMMeshDistribute(DM, const char[], DM*); 575 extern PetscErrorCode DMMeshDistributeByFace(DM, const char[], DM*); 576 extern PetscErrorCode DMMeshGenerate(DM, PetscBool , DM *); 577 extern PetscErrorCode DMMeshRefine(DM, double, PetscBool , DM*); 578 extern PetscErrorCode DMMeshLoad(PetscViewer, DM); 579 extern PetscErrorCode DMMeshGetMaximumDegree(DM, PetscInt *); 580 581 extern PetscErrorCode DMMeshGetLabelSize(DM, const char[], PetscInt *); 582 extern PetscErrorCode DMMeshGetLabelIds(DM, const char[], PetscInt *); 583 extern PetscErrorCode DMMeshGetStratumSize(DM, const char [], PetscInt, PetscInt *); 584 extern PetscErrorCode DMMeshGetStratum(DM, const char [], PetscInt, PetscInt *); 585 586 extern PetscErrorCode DMCartesianCreate(MPI_Comm, DM *); 587 extern PetscErrorCode DMMeshCartesianGetMesh(DM, ALE::Obj<ALE::CartesianMesh>&); 588 extern PetscErrorCode DMMeshCartesianSetMesh(DM, const ALE::Obj<ALE::CartesianMesh>&); 589 590 extern PetscErrorCode DMMeshGetCoordinates(DM, PetscBool , PetscInt *, PetscInt *, PetscReal *[]); 591 extern PetscErrorCode DMMeshGetElements(DM, PetscBool , PetscInt *, PetscInt *, PetscInt *[]); 592 extern PetscErrorCode DMMeshGetCone(DM, PetscInt, PetscInt *, PetscInt *[]); 593 594 extern PetscErrorCode restrictVector(Vec, Vec, InsertMode); 595 extern PetscErrorCode assembleVectorComplete(Vec, Vec, InsertMode); 596 extern PetscErrorCode assembleVector(Vec, PetscInt, PetscScalar [], InsertMode); 597 extern PetscErrorCode updateOperator(Mat, const ALE::Obj<PETSC_MESH_TYPE>&, const ALE::Obj<PETSC_MESH_TYPE::real_section_type>&, const ALE::Obj<PETSC_MESH_TYPE::order_type>&, const PETSC_MESH_TYPE::point_type&, PetscScalar [], InsertMode); 598 extern PetscErrorCode updateOperatorGeneral(Mat, const ALE::Obj<PETSC_MESH_TYPE>&, const ALE::Obj<PETSC_MESH_TYPE::real_section_type>&, const ALE::Obj<PETSC_MESH_TYPE::order_type>&, const PETSC_MESH_TYPE::point_type&, const ALE::Obj<PETSC_MESH_TYPE>&, const ALE::Obj<PETSC_MESH_TYPE::real_section_type>&, const ALE::Obj<PETSC_MESH_TYPE::order_type>&, const PETSC_MESH_TYPE::point_type&, PetscScalar [], InsertMode); 599 600 /*S 601 SectionReal - Abstract PETSc object that manages distributed field data over a topology (Sieve). 602 603 Level: beginner 604 605 Concepts: distributed mesh, field 606 607 .seealso: SectionRealCreate(), SectionRealDestroy(), Mesh, DMMeshCreate() 608 S*/ 609 typedef struct _p_SectionReal* SectionReal; 610 611 /* Logging support */ 612 extern PetscClassId SECTIONREAL_CLASSID; 613 614 extern PetscErrorCode SectionRealCreate(MPI_Comm,SectionReal*); 615 extern PetscErrorCode SectionRealDestroy(SectionReal); 616 extern PetscErrorCode SectionRealView(SectionReal,PetscViewer); 617 extern PetscErrorCode SectionRealDuplicate(SectionReal,SectionReal*); 618 619 extern PetscErrorCode SectionRealGetSection(SectionReal,ALE::Obj<PETSC_MESH_TYPE::real_section_type>&); 620 extern PetscErrorCode SectionRealSetSection(SectionReal,const ALE::Obj<PETSC_MESH_TYPE::real_section_type>&); 621 extern PetscErrorCode SectionRealGetBundle(SectionReal,ALE::Obj<PETSC_MESH_TYPE>&); 622 extern PetscErrorCode SectionRealSetBundle(SectionReal,const ALE::Obj<PETSC_MESH_TYPE>&); 623 624 extern PetscErrorCode SectionRealDistribute(SectionReal, DM, SectionReal *); 625 extern PetscErrorCode SectionRealRestrict(SectionReal, PetscInt, PetscScalar *[]); 626 extern PetscErrorCode SectionRealUpdate(SectionReal, PetscInt, const PetscScalar [], InsertMode); 627 extern PetscErrorCode SectionRealZero(SectionReal); 628 extern PetscErrorCode SectionRealCreateLocalVector(SectionReal, Vec*); 629 extern PetscErrorCode SectionRealAddSpace(SectionReal); 630 extern PetscErrorCode SectionRealGetFibration(SectionReal, const PetscInt, SectionReal *); 631 extern PetscErrorCode SectionRealToVec(SectionReal, DM, ScatterMode, Vec); 632 extern PetscErrorCode SectionRealToVec(SectionReal, VecScatter, ScatterMode, Vec); 633 extern PetscErrorCode SectionRealNorm(SectionReal, DM, NormType, PetscReal *); 634 extern PetscErrorCode SectionRealAXPY(SectionReal, DM, PetscScalar, SectionReal); 635 extern PetscErrorCode SectionRealComplete(SectionReal); 636 extern PetscErrorCode SectionRealSet(SectionReal, PetscReal); 637 extern PetscErrorCode SectionRealGetFiberDimension(SectionReal, PetscInt, PetscInt*); 638 extern PetscErrorCode SectionRealSetFiberDimension(SectionReal, PetscInt, const PetscInt); 639 extern PetscErrorCode SectionRealSetFiberDimensionField(SectionReal, PetscInt, const PetscInt, const PetscInt); 640 extern PetscErrorCode SectionRealGetSize(SectionReal, PetscInt *); 641 extern PetscErrorCode SectionRealAllocate(SectionReal); 642 extern PetscErrorCode SectionRealClear(SectionReal); 643 644 extern PetscErrorCode SectionRealRestrictClosure(SectionReal, DM, PetscInt, PetscInt, PetscScalar []); 645 extern PetscErrorCode SectionRealRestrictClosure(SectionReal, DM, PetscInt, const PetscScalar *[]); 646 extern PetscErrorCode SectionRealUpdateClosure(SectionReal, DM, PetscInt, PetscScalar [], InsertMode); 647 648 extern PetscErrorCode DMMeshHasSectionReal(DM, const char [], PetscBool *); 649 extern PetscErrorCode DMMeshGetSectionReal(DM, const char [], SectionReal *); 650 extern PetscErrorCode DMMeshSetSectionReal(DM, SectionReal); 651 extern PetscErrorCode DMMeshCreateMatrix(DM, SectionReal, MatType, Mat *); 652 extern PetscErrorCode DMMeshCreateVector(DM, SectionReal, Vec *); 653 extern PetscErrorCode DMMeshCreateGlobalScatter(DM, SectionReal, VecScatter *); 654 extern PetscErrorCode assembleVector(Vec, DM, SectionReal, PetscInt, PetscScalar [], InsertMode); 655 extern PetscErrorCode assembleMatrix(Mat, DM, SectionReal, PetscInt, PetscScalar [], InsertMode); 656 657 extern PetscErrorCode DMMeshCreateGlobalRealVector(DM, SectionReal, Vec *); 658 extern PetscErrorCode DMMeshGetGlobalScatter(DM,VecScatter *); 659 extern PetscErrorCode DMMeshCreateGlobalScatter(DM,SectionReal,VecScatter *); 660 extern PetscErrorCode DMMeshGetLocalFunction(DM, PetscErrorCode (**)(DM, SectionReal, SectionReal, void*)); 661 extern PetscErrorCode DMMeshSetLocalFunction(DM, PetscErrorCode (*)(DM, SectionReal, SectionReal, void*)); 662 extern PetscErrorCode DMMeshGetLocalJacobian(DM, PetscErrorCode (**)(DM, SectionReal, Mat, void*)); 663 extern PetscErrorCode DMMeshSetLocalJacobian(DM, PetscErrorCode (*)(DM, SectionReal, Mat, void*)); 664 extern PetscErrorCode DMMeshFormFunction(DM, SectionReal, SectionReal, void*); 665 extern PetscErrorCode DMMeshFormJacobian(DM, SectionReal, Mat, void*); 666 extern PetscErrorCode DMMeshInterpolatePoints(DM, SectionReal, int, double *, double **); 667 668 /*S 669 SectionInt - Abstract PETSc object that manages distributed field data over a topology (Sieve). 670 671 Level: beginner 672 673 Concepts: distributed mesh, field 674 675 .seealso: SectionIntCreate(), SectionIntDestroy(), DM, DMMeshCreate() 676 S*/ 677 typedef struct _p_SectionInt* SectionInt; 678 679 /* Logging support */ 680 extern PetscClassId SECTIONINT_CLASSID; 681 682 extern PetscErrorCode SectionIntCreate(MPI_Comm,SectionInt*); 683 extern PetscErrorCode SectionIntDestroy(SectionInt); 684 extern PetscErrorCode SectionIntView(SectionInt,PetscViewer); 685 686 extern PetscErrorCode SectionIntGetSection(SectionInt,ALE::Obj<PETSC_MESH_TYPE::int_section_type>&); 687 extern PetscErrorCode SectionIntSetSection(SectionInt,const ALE::Obj<PETSC_MESH_TYPE::int_section_type>&); 688 extern PetscErrorCode SectionIntGetBundle(SectionInt,ALE::Obj<PETSC_MESH_TYPE>&); 689 extern PetscErrorCode SectionIntSetBundle(SectionInt,const ALE::Obj<PETSC_MESH_TYPE>&); 690 691 extern PetscErrorCode SectionIntDistribute(SectionInt, DM, SectionInt *); 692 extern PetscErrorCode SectionIntRestrict(SectionInt, PetscInt, PetscInt *[]); 693 extern PetscErrorCode SectionIntUpdate(SectionInt, PetscInt, const PetscInt [], InsertMode); 694 extern PetscErrorCode SectionIntZero(SectionInt); 695 extern PetscErrorCode SectionIntComplete(SectionInt); 696 extern PetscErrorCode SectionIntGetFiberDimension(SectionInt, PetscInt, PetscInt*); 697 extern PetscErrorCode SectionIntSetFiberDimension(SectionInt, PetscInt, const PetscInt); 698 extern PetscErrorCode SectionIntSetFiberDimensionField(SectionInt, PetscInt, const PetscInt, const PetscInt); 699 extern PetscErrorCode SectionIntGetSize(SectionInt, PetscInt *); 700 extern PetscErrorCode SectionIntAllocate(SectionInt); 701 extern PetscErrorCode SectionIntClear(SectionInt); 702 703 extern PetscErrorCode SectionIntAddSpace(SectionInt); 704 extern PetscErrorCode SectionIntGetFibration(SectionInt, const PetscInt, SectionInt *); 705 extern PetscErrorCode SectionIntSet(SectionInt, PetscInt); 706 707 extern PetscErrorCode SectionIntRestrictClosure(SectionInt, DM, PetscInt, PetscInt, PetscInt []); 708 extern PetscErrorCode SectionIntUpdateClosure(SectionInt, DM, PetscInt, PetscInt [], InsertMode); 709 710 extern PetscErrorCode DMMeshHasSectionInt(DM, const char [], PetscBool *); 711 extern PetscErrorCode DMMeshGetSectionInt(DM, const char [], SectionInt *); 712 extern PetscErrorCode DMMeshSetSectionInt(DM, SectionInt); 713 714 typedef PetscErrorCode (*DMMeshLocalFunction1)(DM,SectionReal,SectionReal,void*); 715 typedef PetscErrorCode (*DMMeshLocalJacobian1)(DM,SectionReal,Mat,void*); 716 717 /* Misc Mesh functions*/ 718 extern PetscErrorCode DMMeshSetMaxDof(DM, PetscInt); 719 720 /* Helper functions for simple distributions */ 721 extern PetscErrorCode DMMeshGetVertexMatrix(DM, MatType, Mat *); 722 extern PetscErrorCode DMMeshGetVertexSectionReal(DM, const char[], PetscInt, SectionReal *); 723 PetscPolymorphicSubroutine(DMMeshGetVertexSectionReal,(DM dm, PetscInt fiberDim, SectionReal *section),(dm,"default",fiberDim,section)) 724 extern PetscErrorCode DMMeshGetVertexSectionInt(DM, const char[], PetscInt, SectionInt *); 725 PetscPolymorphicSubroutine(DMMeshGetVertexSectionInt,(DM dm, PetscInt fiberDim, SectionInt *section),(dm,"default",fiberDim,section)) 726 extern PetscErrorCode DMMeshGetCellMatrix(DM, MatType, Mat *); 727 extern PetscErrorCode DMMeshGetCellSectionReal(DM, const char[], PetscInt, SectionReal *); 728 PetscPolymorphicSubroutine(DMMeshGetCellSectionReal,(DM dm, PetscInt fiberDim, SectionReal *section),(dm,"default",fiberDim,section)) 729 extern PetscErrorCode DMMeshGetCellSectionInt(DM, const char[], PetscInt, SectionInt *); 730 PetscPolymorphicSubroutine(DMMeshGetCellSectionInt,(DM dm, PetscInt fiberDim, SectionInt *section),(dm,"default",fiberDim,section)) 731 732 /* Support for various mesh formats */ 733 extern PetscErrorCode DMMeshCreateExodus(MPI_Comm, const char [], DM *); 734 extern PetscErrorCode DMMeshExodusGetInfo(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt *); 735 736 extern PetscErrorCode DMWriteVTKHeader(PetscViewer); 737 extern PetscErrorCode DMWriteVTKVertices(DM, PetscViewer); 738 extern PetscErrorCode DMWriteVTKElements(DM, PetscViewer); 739 extern PetscErrorCode DMWritePCICEVertices(DM, PetscViewer); 740 extern PetscErrorCode DMWritePCICEElements(DM, PetscViewer); 741 extern PetscErrorCode DMWritePyLithVertices(DM, PetscViewer); 742 extern PetscErrorCode DMWritePyLithElements(DM, SectionReal, PetscViewer); 743 extern PetscErrorCode DMWritePyLithVerticesLocal(DM, PetscViewer); 744 extern PetscErrorCode WDMritePyLithElementsLocal(DM, SectionReal, PetscViewer); 745 746 typedef struct { 747 int numQuadPoints, numBasisFuncs; 748 const double *quadPoints, *quadWeights, *basis, *basisDer; 749 } PetscQuadrature; 750 751 #endif /* Mesh section */ 752 753 PETSC_EXTERN_CXX_END 754 #endif 755