Actual source code: petscmat.h
1: /* $Id: petscmat.h,v 1.228 2001/09/07 20:09:08 bsmith Exp $ */
2: /*
3: Include file for the matrix component of PETSc
4: */
5: #ifndef __PETSCMAT_H
7: #include petscvec.h
9: /*S
10: Mat - Abstract PETSc matrix object
12: Level: beginner
14: Concepts: matrix; linear operator
16: .seealso: MatCreate(), MatType, MatSetType()
17: S*/
18: typedef struct _p_Mat* Mat;
20: /*E
21: MatType - String with the name of a PETSc matrix or the creation function
22: with an optional dynamic library name, for example
23: http://www.mcs.anl.gov/petsc/lib.a:mymatcreate()
25: Level: beginner
27: .seealso: MatSetType(), Mat
28: E*/
29: #define MATSAME "same"
30: #define MATSEQMAIJ "seqmaij"
31: #define MATMPIMAIJ "mpimaij"
32: #define MATIS "is"
33: #define MATMPIROWBS "mpirowbs"
34: #define MATSEQDENSE "seqdense"
35: #define MATSEQAIJ "seqaij"
36: #define MATMPIAIJ "mpiaij"
37: #define MATSHELL "shell"
38: #define MATSEQBDIAG "seqbdiag"
39: #define MATMPIBDIAG "mpibdiag"
40: #define MATMPIDENSE "mpidense"
41: #define MATSEQBAIJ "seqbaij"
42: #define MATMPIBAIJ "mpibaij"
43: #define MATMPIADJ "mpiadj"
44: #define MATSEQSBAIJ "seqsbaij"
45: #define MATMPISBAIJ "mpisbaij"
46: #define MATDAAD "daad"
47: #define MATMFFD "mffd"
48: #define MATESI "esi"
49: #define MATPETSCESI "petscesi"
50: typedef char* MatType;
52: #define MAT_SER_SEQAIJ_BINARY "seqaij_binary"
53: #define MAT_SER_MPIAIJ_BINARY "mpiaij_binary"
54: typedef char *MatSerializeType;
56: /* Logging support */
57: #define MAT_FILE_COOKIE 1211216 /* used to indicate matrices in binary files */
58: extern int MAT_COOKIE;
59: extern int MATSNESMFCTX_COOKIE;
60: extern int MAT_FDCOLORING_COOKIE;
61: extern int MAT_PARTITIONING_COOKIE;
62: extern int MAT_NULLSPACE_COOKIE;
63: extern int MAT_Mult, MAT_MultMatrixFree, MAT_MultMultiple, MAT_MultConstrained, MAT_MultAdd, MAT_MultTranspose;
64: extern int MAT_MultTransposeConstrained, MAT_MultTransposeAdd, MAT_Solve, MAT_SolveMultiple, MAT_SolveAdd, MAT_SolveTranspose;
65: extern int MAT_SolveTransposeAdd, MAT_Relax, MAT_ForwardSolve, MAT_BackwardSolve, MAT_LUFactor, MAT_LUFactorSymbolic;
66: extern int MAT_LUFactorNumeric, MAT_CholeskyFactor, MAT_CholeskyFactorSymbolic, MAT_CholeskyFactorNumeric, MAT_ILUFactor;
67: extern int MAT_ILUFactorSymbolic, MAT_ICCFactorSymbolic, MAT_Copy, MAT_Convert, MAT_Scale, MAT_AssemblyBegin;
68: extern int MAT_AssemblyEnd, MAT_SetValues, MAT_GetValues, MAT_GetRow, MAT_GetSubMatrices, MAT_GetColoring, MAT_GetOrdering;
69: extern int MAT_IncreaseOverlap, MAT_Partitioning, MAT_ZeroEntries, MAT_Load, MAT_View, MAT_AXPY, MAT_FDColoringCreate;
70: extern int MAT_FDColoringApply, MAT_Transpose;
72: EXTERN int MatInitializePackage(char *);
74: EXTERN int MatCreate(MPI_Comm,int,int,int,int,Mat*);
75: EXTERN int MatSetType(Mat,MatType);
76: EXTERN int MatSetFromOptions(Mat);
77: EXTERN int MatSetUpPreallocation(Mat);
78: EXTERN int MatRegisterAll(char*);
79: EXTERN int MatRegister(char*,char*,char*,int(*)(Mat));
80: EXTERN int MatSerializeRegister(const char [], const char [], const char [], int (*)(MPI_Comm, Mat *, PetscViewer, PetscTruth));
81: #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
82: #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,0)
83: #define MatSerializeRegisterDynamic(a,b,c,d) MatSerializeRegister(a,b,c,0)
84: #else
85: #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,d)
86: #define MatSerializeRegisterDynamic(a,b,c,d) MatSerializeRegister(a,b,c,d)
87: #endif
88: extern PetscTruth MatRegisterAllCalled;
89: extern PetscFList MatList;
91: EXTERN PetscFList MatSerializeList;
92: EXTERN int MatSerializeRegisterAll(const char []);
93: EXTERN int MatSerializeRegisterDestroy(void);
94: EXTERN int MatSerializeRegisterAllCalled;
95: EXTERN int MatSerialize(MPI_Comm, Mat *, PetscViewer, PetscTruth);
96: EXTERN int MatSetSerializeType(Mat, MatSerializeType);
98: EXTERN int MatCreate(MPI_Comm,int,int,int,int,Mat*);
99: EXTERN int MatCreateSeqDense(MPI_Comm,int,int,PetscScalar*,Mat*);
100: EXTERN int MatCreateMPIDense(MPI_Comm,int,int,int,int,PetscScalar*,Mat*);
101: EXTERN int MatCreateSeqAIJ(MPI_Comm,int,int,int,int*,Mat*);
102: EXTERN int MatCreateMPIAIJ(MPI_Comm,int,int,int,int,int,int*,int,int*,Mat*);
103: EXTERN int MatCreateMPIRowbs(MPI_Comm,int,int,int,int*,Mat*);
104: EXTERN int MatCreateSeqBDiag(MPI_Comm,int,int,int,int,int*,PetscScalar**,Mat*);
105: EXTERN int MatCreateMPIBDiag(MPI_Comm,int,int,int,int,int,int*,PetscScalar**,Mat*);
106: EXTERN int MatCreateSeqBAIJ(MPI_Comm,int,int,int,int,int*,Mat*);
107: EXTERN int MatCreateMPIBAIJ(MPI_Comm,int,int,int,int,int,int,int*,int,int*,Mat*);
108: EXTERN int MatCreateMPIAdj(MPI_Comm,int,int,int*,int*,int *,Mat*);
109: EXTERN int MatCreateSeqSBAIJ(MPI_Comm,int,int,int,int,int*,Mat*);
110: EXTERN int MatCreateMPISBAIJ(MPI_Comm,int,int,int,int,int,int,int*,int,int*,Mat*);
111: EXTERN int MatCreateShell(MPI_Comm,int,int,int,int,void *,Mat*);
112: EXTERN int MatCreateAdic(MPI_Comm,int,int,int,int,int,void (*)(void),Mat*);
113: EXTERN int MatDestroy(Mat);
115: EXTERN int MatPrintHelp(Mat);
116: EXTERN int MatGetPetscMaps(Mat,PetscMap*,PetscMap*);
118: /* ------------------------------------------------------------*/
119: EXTERN int MatSetValues(Mat,int,int*,int,int*,PetscScalar*,InsertMode);
120: EXTERN int MatSetValuesBlocked(Mat,int,int*,int,int*,PetscScalar*,InsertMode);
122: /*S
123: MatStencil - Data structure (C struct) for storing information about a single row or
124: column of a matrix as index on an associated grid.
126: Level: beginner
128: Concepts: matrix; linear operator
130: .seealso: MatSetValuesStencil(), MatSetStencil()
131: S*/
132: typedef struct {
133: int k,j,i,c;
134: } MatStencil;
136: EXTERN int MatSetValuesStencil(Mat,int,MatStencil*,int,MatStencil*,PetscScalar*,InsertMode);
137: EXTERN int MatSetValuesBlockedStencil(Mat,int,MatStencil*,int,MatStencil*,PetscScalar*,InsertMode);
138: EXTERN int MatSetStencil(Mat,int,int*,int*,int);
140: EXTERN int MatSetColoring(Mat,ISColoring);
141: EXTERN int MatSetValuesAdic(Mat,void*);
142: EXTERN int MatSetValuesAdifor(Mat,int,void*);
144: /*E
145: MatAssemblyType - Indicates if the matrix is now to be used, or if you plan
146: to continue to add values to it
148: Level: beginner
150: .seealso: MatAssemblyBegin(), MatAssemblyEnd()
151: E*/
152: typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType;
153: EXTERN int MatAssemblyBegin(Mat,MatAssemblyType);
154: EXTERN int MatAssemblyEnd(Mat,MatAssemblyType);
155: EXTERN int MatAssembled(Mat,PetscTruth*);
157: #define MatSetValue(v,i,j,va,mode)
158: 0; {int _ierr,_row = i,_col = j; PetscScalar _va = va;
159: _MatSetValues(v,1,&_row,1,&_col,&_va,mode);CHKERRQ(_ierr);
160: }
161: #define MatGetValue(v,i,j,va)
162: 0; {int _ierr,_row = i,_col = j;
163: _MatGetValues(v,1,&_row,1,&_col,&va);CHKERRQ(_ierr);
164: }
165: #define MatSetValueLocal(v,i,j,va,mode)
166: 0; {int _ierr,_row = i,_col = j; PetscScalar _va = va;
167: _MatSetValuesLocal(v,1,&_row,1,&_col,&_va,mode);CHKERRQ(_ierr);
168: }
169: /*E
170: MatOption - Options that may be set for a matrix and its behavior or storage
172: Level: beginner
174: Any additions/changes here MUST also be made in include/finclude/petscmat.h
176: .seealso: MatSetOption()
177: E*/
178: typedef enum {MAT_ROW_ORIENTED=1,MAT_COLUMN_ORIENTED=2,MAT_ROWS_SORTED=4,
179: MAT_COLUMNS_SORTED=8,MAT_NO_NEW_NONZERO_LOCATIONS=16,
180: MAT_YES_NEW_NONZERO_LOCATIONS=32,MAT_SYMMETRIC=64,
181: MAT_STRUCTURALLY_SYMMETRIC=65,MAT_NO_NEW_DIAGONALS=66,
182: MAT_YES_NEW_DIAGONALS=67,MAT_INODE_LIMIT_1=68,MAT_INODE_LIMIT_2=69,
183: MAT_INODE_LIMIT_3=70,MAT_INODE_LIMIT_4=71,MAT_INODE_LIMIT_5=72,
184: MAT_IGNORE_OFF_PROC_ENTRIES=73,MAT_ROWS_UNSORTED=74,
185: MAT_COLUMNS_UNSORTED=75,MAT_NEW_NONZERO_LOCATION_ERR=76,
186: MAT_NEW_NONZERO_ALLOCATION_ERR=77,MAT_USE_HASH_TABLE=78,
187: MAT_KEEP_ZEROED_ROWS=79,MAT_IGNORE_ZERO_ENTRIES=80,MAT_USE_INODES=81,
188: MAT_DO_NOT_USE_INODES=82,MAT_USE_SINGLE_PRECISION_SOLVES=83} MatOption;
189: EXTERN int MatSetOption(Mat,MatOption);
190: EXTERN int MatGetType(Mat,MatType*);
192: EXTERN int MatGetValues(Mat,int,int*,int,int*,PetscScalar*);
193: EXTERN int MatGetRow(Mat,int,int *,int **,PetscScalar**);
194: EXTERN int MatRestoreRow(Mat,int,int *,int **,PetscScalar**);
195: EXTERN int MatGetColumn(Mat,int,int *,int **,PetscScalar**);
196: EXTERN int MatRestoreColumn(Mat,int,int *,int **,PetscScalar**);
197: EXTERN int MatGetColumnVector(Mat,Vec,int);
198: EXTERN int MatGetArray(Mat,PetscScalar **);
199: EXTERN int MatRestoreArray(Mat,PetscScalar **);
200: EXTERN int MatGetBlockSize(Mat,int *);
202: EXTERN int MatMult(Mat,Vec,Vec);
203: EXTERN int MatMultAdd(Mat,Vec,Vec,Vec);
204: EXTERN int MatMultTranspose(Mat,Vec,Vec);
205: EXTERN int MatMultTransposeAdd(Mat,Vec,Vec,Vec);
206: EXTERN int MatMultConstrained(Mat,Vec,Vec);
207: EXTERN int MatMultTransposeConstrained(Mat,Vec,Vec);
209: /*E
210: MatDuplicateOption - Indicates if a duplicated sparse matrix should have
211: its numerical values copied over or just its nonzero structure.
213: Level: beginner
215: Any additions/changes here MUST also be made in include/finclude/petscmat.h
217: .seealso: MatDuplicate()
218: E*/
219: typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES} MatDuplicateOption;
221: EXTERN int MatConvertRegister(char*,char*,char*,int (*)(Mat,MatType,Mat*));
222: #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
223: #define MatConvertRegisterDynamic(a,b,c,d) MatConvertRegister(a,b,c,0)
224: #else
225: #define MatConvertRegisterDynamic(a,b,c,d) MatConvertRegister(a,b,c,d)
226: #endif
227: EXTERN int MatConvertRegisterAll(char*);
228: EXTERN int MatConvertRegisterDestroy(void);
229: extern PetscTruth MatConvertRegisterAllCalled;
230: extern PetscFList MatConvertList;
231: EXTERN int MatConvert(Mat,MatType,Mat*);
232: EXTERN int MatDuplicate(Mat,MatDuplicateOption,Mat*);
234: /*E
235: MatStructure - Indicates if the matrix has the same nonzero structure
237: Level: beginner
239: Any additions/changes here MUST also be made in include/finclude/petscmat.h
241: .seealso: MatCopy(), SLESSetOperators(), PCSetOperators()
242: E*/
243: typedef enum {SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER} MatStructure;
245: EXTERN int MatCopy(Mat,Mat,MatStructure);
246: EXTERN int MatView(Mat,PetscViewer);
248: EXTERN int MatLoadRegister(char*,char*,char*,int (*)(PetscViewer,MatType,Mat*));
249: #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
250: #define MatLoadRegisterDynamic(a,b,c,d) MatLoadRegister(a,b,c,0)
251: #else
252: #define MatLoadRegisterDynamic(a,b,c,d) MatLoadRegister(a,b,c,d)
253: #endif
254: EXTERN int MatLoadRegisterAll(char*);
255: EXTERN int MatLoadRegisterDestroy(void);
256: extern PetscTruth MatLoadRegisterAllCalled;
257: extern PetscFList MatLoadList;
258: EXTERN int MatLoad(PetscViewer,MatType,Mat*);
260: EXTERN int MatGetRowIJ(Mat,int,PetscTruth,int*,int **,int **,PetscTruth *);
261: EXTERN int MatRestoreRowIJ(Mat,int,PetscTruth,int *,int **,int **,PetscTruth *);
262: EXTERN int MatGetColumnIJ(Mat,int,PetscTruth,int*,int **,int **,PetscTruth *);
263: EXTERN int MatRestoreColumnIJ(Mat,int,PetscTruth,int *,int **,int **,PetscTruth *);
265: /*S
266: MatInfo - Context of matrix information, used with MatGetInfo()
268: In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE
270: Level: intermediate
272: Concepts: matrix^nonzero information
274: .seealso: MatGetInfo(), MatInfoType
275: S*/
276: typedef struct {
277: PetscLogDouble rows_global,columns_global; /* number of global rows and columns */
278: PetscLogDouble rows_local,columns_local; /* number of local rows and columns */
279: PetscLogDouble block_size; /* block size */
280: PetscLogDouble nz_allocated,nz_used,nz_unneeded; /* number of nonzeros */
281: PetscLogDouble memory; /* memory allocated */
282: PetscLogDouble assemblies; /* number of matrix assemblies called */
283: PetscLogDouble mallocs; /* number of mallocs during MatSetValues() */
284: PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
285: PetscLogDouble factor_mallocs; /* number of mallocs during factorization */
286: } MatInfo;
288: /*E
289: MatInfoType - Indicates if you want information about the local part of the matrix,
290: the entire parallel matrix or the maximum over all the local parts.
292: Level: beginner
294: Any additions/changes here MUST also be made in include/finclude/petscmat.h
296: .seealso: MatGetInfo(), MatInfo
297: E*/
298: typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
299: EXTERN int MatGetInfo(Mat,MatInfoType,MatInfo*);
300: EXTERN int MatValid(Mat,PetscTruth*);
301: EXTERN int MatGetDiagonal(Mat,Vec);
302: EXTERN int MatGetRowMax(Mat,Vec);
303: EXTERN int MatTranspose(Mat,Mat*);
304: EXTERN int MatPermute(Mat,IS,IS,Mat *);
305: EXTERN int MatPermuteSparsify(Mat,int,PetscReal,PetscReal,IS,IS,Mat *);
306: EXTERN int MatDiagonalScale(Mat,Vec,Vec);
307: EXTERN int MatDiagonalSet(Mat,Vec,InsertMode);
308: EXTERN int MatEqual(Mat,Mat,PetscTruth*);
310: EXTERN int MatNorm(Mat,NormType,PetscReal *);
311: EXTERN int MatZeroEntries(Mat);
312: EXTERN int MatZeroRows(Mat,IS,PetscScalar*);
313: EXTERN int MatZeroColumns(Mat,IS,PetscScalar*);
315: EXTERN int MatUseScaledForm(Mat,PetscTruth);
316: EXTERN int MatScaleSystem(Mat,Vec,Vec);
317: EXTERN int MatUnScaleSystem(Mat,Vec,Vec);
319: EXTERN int MatGetSize(Mat,int*,int*);
320: EXTERN int MatGetLocalSize(Mat,int*,int*);
321: EXTERN int MatGetOwnershipRange(Mat,int*,int*);
323: /*E
324: MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices()
325: or MatGetSubMatrix() are to be reused to store the new matrix values.
327: Level: beginner
329: Any additions/changes here MUST also be made in include/finclude/petscmat.h
331: .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices()
332: E*/
333: typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX} MatReuse;
334: EXTERN int MatGetSubMatrices(Mat,int,IS *,IS *,MatReuse,Mat **);
335: EXTERN int MatDestroyMatrices(int,Mat **);
336: EXTERN int MatGetSubMatrix(Mat,IS,IS,int,MatReuse,Mat *);
338: EXTERN int MatIncreaseOverlap(Mat,int,IS *,int);
340: EXTERN int MatAXPY(PetscScalar *,Mat,Mat,MatStructure);
341: EXTERN int MatAYPX(PetscScalar *,Mat,Mat);
342: EXTERN int MatCompress(Mat);
344: EXTERN int MatScale(PetscScalar *,Mat);
345: EXTERN int MatShift(PetscScalar *,Mat);
347: EXTERN int MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping);
348: EXTERN int MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping);
349: EXTERN int MatZeroRowsLocal(Mat,IS,PetscScalar*);
350: EXTERN int MatSetValuesLocal(Mat,int,int*,int,int*,PetscScalar*,InsertMode);
351: EXTERN int MatSetValuesBlockedLocal(Mat,int,int*,int,int*,PetscScalar*,InsertMode);
353: EXTERN int MatSetStashInitialSize(Mat,int,int);
355: EXTERN int MatInterpolateAdd(Mat,Vec,Vec,Vec);
356: EXTERN int MatInterpolate(Mat,Vec,Vec);
357: EXTERN int MatRestrict(Mat,Vec,Vec);
359: /*
360: These three (or four) macros MUST be used together. The third one closes the open { of the first one
361: */
362: #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0;
363: {
364: int _4_ierr,__tmp = (nrows),__ctmp = (ncols),__rstart,__start,__end;
365: _4_PetscMalloc(2*__tmp*sizeof(int),&dnz);CHKERRQ(_4_ierr);onz = dnz + __tmp;
366: _4_PetscMemzero(dnz,2*__tmp*sizeof(int));CHKERRQ(_4_ierr);
367: _4_MPI_Scan(&__ctmp,&__end,1,MPI_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;
368: _4_MPI_Scan(&__tmp,&__rstart,1,MPI_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __tmp;
370: #define MatPreallocateSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;
371: {
372: int __l;
373: _4_ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);
374: _4_ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);
375: for (__l=0;__l<nrows;__l++) {
376: _4_MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);
377: }
378: }
379:
380: #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;
381: { int __i;
382: for (__i=0; __i<nc; __i++) {
383: if (cols[__i] < __start || cols[__i] >= __end) onz[row - __rstart]++;
384: }
385: dnz[row - __rstart] = nc - onz[row - __rstart];
386: }
388: #define MatPreallocateFinalize(dnz,onz) 0;_4_PetscFree(dnz);CHKERRQ(_4_ierr);}
390: /* Routines unique to particular data structures */
391: EXTERN int MatShellGetContext(Mat,void **);
393: EXTERN int MatBDiagGetData(Mat,int*,int*,int**,int**,PetscScalar***);
394: EXTERN int MatSeqAIJSetColumnIndices(Mat,int *);
395: EXTERN int MatSeqBAIJSetColumnIndices(Mat,int *);
396: EXTERN int MatCreateSeqAIJWithArrays(MPI_Comm,int,int,int*,int*,PetscScalar *,Mat*);
398: EXTERN int MatSeqBAIJSetPreallocation(Mat,int,int,int*);
399: EXTERN int MatSeqSBAIJSetPreallocation(Mat,int,int,int*);
400: EXTERN int MatSeqAIJSetPreallocation(Mat,int,int*);
401: EXTERN int MatSeqDensePreallocation(Mat,PetscScalar*);
402: EXTERN int MatSeqBDiagSetPreallocation(Mat,int,int,int*,PetscScalar**);
403: EXTERN int MatSeqDenseSetPreallocation(Mat,PetscScalar*);
405: EXTERN int MatMPIBAIJSetPreallocation(Mat,int,int,int*,int,int*);
406: EXTERN int MatMPISBAIJSetPreallocation(Mat,int,int,int*,int,int*);
407: EXTERN int MatMPIAIJSetPreallocation(Mat,int,int*,int,int*);
408: EXTERN int MatMPIDensePreallocation(Mat,PetscScalar*);
409: EXTERN int MatMPIBDiagSetPreallocation(Mat,int,int,int*,PetscScalar**);
410: EXTERN int MatMPIAdjSetPreallocation(Mat,int*,int*,int*);
411: EXTERN int MatMPIDenseSetPreallocation(Mat,PetscScalar*);
412: EXTERN int MatMPIRowbsSetPreallocation(Mat,int,int*);
413: EXTERN int MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,int**);
414: EXTERN int MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,int**);
415: EXTERN int MatAdicSetLocalFunction(Mat,void (*)(void));
417: EXTERN int MatStoreValues(Mat);
418: EXTERN int MatRetrieveValues(Mat);
420: EXTERN int MatDAADSetCtx(Mat,void*);
422: /*
423: These routines are not usually accessed directly, rather solving is
424: done through the SLES, KSP and PC interfaces.
425: */
427: /*E
428: MatOrderingType - String with the name of a PETSc matrix ordering or the creation function
429: with an optional dynamic library name, for example
430: http://www.mcs.anl.gov/petsc/lib.a:orderingcreate()
432: Level: beginner
434: .seealso: MatGetOrdering()
435: E*/
436: typedef char* MatOrderingType;
437: #define MATORDERING_NATURAL "natural"
438: #define MATORDERING_ND "nd"
439: #define MATORDERING_1WD "1wd"
440: #define MATORDERING_RCM "rcm"
441: #define MATORDERING_QMD "qmd"
442: #define MATORDERING_ROWLENGTH "rowlength"
443: #define MATORDERING_DSC_ND "dsc_nd"
444: #define MATORDERING_DSC_MMD "dsc_mmd"
445: #define MATORDERING_DSC_MDF "dsc_mdf"
446: #define MATORDERING_CONSTRAINED "constrained"
447: #define MATORDERING_IDENTITY "identity"
448: #define MATORDERING_REVERSE "reverse"
450: EXTERN int MatGetOrdering(Mat,MatOrderingType,IS*,IS*);
451: EXTERN int MatOrderingRegister(char*,char*,char*,int(*)(Mat,MatOrderingType,IS*,IS*));
452: #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
453: #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,0)
454: #else
455: #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,d)
456: #endif
457: EXTERN int MatOrderingRegisterDestroy(void);
458: EXTERN int MatOrderingRegisterAll(char*);
459: extern PetscTruth MatOrderingRegisterAllCalled;
460: extern PetscFList MatOrderingList;
462: EXTERN int MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
464: EXTERN int MatCholeskyFactor(Mat,IS,PetscReal);
465: EXTERN int MatCholeskyFactorSymbolic(Mat,IS,PetscReal,Mat*);
466: EXTERN int MatCholeskyFactorNumeric(Mat,Mat*);
468: /*S
469: MatILUInfo - Data based into the matrix ILU factorization routines
471: In Fortran these are simply double precision arrays of size MAT_ILUINFO_SIZE
473: Notes: These are not usually directly used by users, instead use the PC type of ILU
474: All entries are double precision.
476: Level: developer
478: .seealso: MatILUFactorSymbolic(), MatILUFactor(), MatLUInfo, MatCholeskyInfo
480: S*/
481: typedef struct {
482: PetscReal levels; /* ILU(levels) */
483: PetscReal fill; /* expected fill; nonzeros in factored matrix/nonzeros in original matrix*/
484: PetscReal diagonal_fill; /* force diagonal to fill in if initially not filled */
485: PetscReal dt; /* drop tolerance */
486: PetscReal dtcol; /* tolerance for pivoting */
487: PetscReal dtcount; /* maximum nonzeros to be allowed per row */
488: PetscReal damping; /* scaling of identity added to matrix to prevent zero pivots */
489: PetscReal damp; /* if is 1.0 and factorization fails, damp until successful */
490: PetscReal zeropivot; /* pivot is called zero if less than this */
491: PetscReal pivotinblocks; /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
492: factorization may be faster if do not pivot */
493: } MatILUInfo;
495: /*S
496: MatLUInfo - Data based into the matrix LU factorization routines
498: In Fortran these are simply double precision arrays of size MAT_LUINFO_SIZE
500: Notes: These are not usually directly used by users, instead use the PC type of LU
501: All entries are double precision.
503: Level: developer
505: .seealso: MatLUFactorSymbolic(), MatILUInfo, MatCholeskyInfo
507: S*/
508: typedef struct {
509: PetscReal fill; /* expected fill; nonzeros in factored matrix/nonzeros in original matrix */
510: PetscReal dtcol; /* tolerance for pivoting; pivot if off_diagonal*dtcol > diagonal */
511: PetscReal damping; /* scaling of identity added to matrix to prevent zero pivots */
512: PetscReal damp; /* if this is 1.0 and factorization fails, damp until successful */
513: PetscReal zeropivot; /* pivot is called zero if less than this */
514: PetscReal pivotinblocks; /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
515: factorization may be faster if do not pivot */
516: } MatLUInfo;
518: /*S
519: MatCholeskyInfo - Data based into the matrix Cholesky factorization routines
521: In Fortran these are simply double precision arrays of size MAT_CHOLESKYINFO_SIZE
523: Notes: These are not usually directly used by users, instead use the PC type of Cholesky
524: All entries are double precision.
526: Level: developer
528: .seealso: MatCholeskyFactorSymbolic(), MatLUInfo, MatILUInfo
530: S*/
531: typedef struct {
532: PetscReal fill; /* expected fill; nonzeros in factored matrix/nonzeros in original matrix */
533: PetscReal damping; /* scaling of identity added to matrix to prevent zero pivots */
534: PetscReal damp; /* if this is 1.0 and factorization fails, damp until successful */
535: PetscReal pivotinblocks; /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
536: factorization may be faster if do not pivot */
537: } MatCholeskyInfo;
539: EXTERN int MatLUFactor(Mat,IS,IS,MatLUInfo*);
540: EXTERN int MatILUFactor(Mat,IS,IS,MatILUInfo*);
541: EXTERN int MatLUFactorSymbolic(Mat,IS,IS,MatLUInfo*,Mat*);
542: EXTERN int MatILUFactorSymbolic(Mat,IS,IS,MatILUInfo*,Mat*);
543: EXTERN int MatICCFactorSymbolic(Mat,IS,PetscReal,int,Mat*);
544: EXTERN int MatICCFactor(Mat,IS,PetscReal,int);
545: EXTERN int MatLUFactorNumeric(Mat,Mat*);
546: EXTERN int MatILUDTFactor(Mat,MatILUInfo*,IS,IS,Mat *);
548: EXTERN int MatSolve(Mat,Vec,Vec);
549: EXTERN int MatForwardSolve(Mat,Vec,Vec);
550: EXTERN int MatBackwardSolve(Mat,Vec,Vec);
551: EXTERN int MatSolveAdd(Mat,Vec,Vec,Vec);
552: EXTERN int MatSolveTranspose(Mat,Vec,Vec);
553: EXTERN int MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
555: EXTERN int MatSetUnfactored(Mat);
557: /* MatSORType may be bitwise ORd together, so do not change the numbers */
558: /*E
559: MatSORType - What type of (S)SOR to perform
561: Level: beginner
563: May be bitwise ORd together
565: Any additions/changes here MUST also be made in include/finclude/petscmat.h
567: .seealso: MatRelax()
568: E*/
569: typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
570: SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
571: SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
572: SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
573: EXTERN int MatRelax(Mat,Vec,PetscReal,MatSORType,PetscReal,int,int,Vec);
575: /*
576: These routines are for efficiently computing Jacobians via finite differences.
577: */
579: /*E
580: MatColoringType - String with the name of a PETSc matrix coloring or the creation function
581: with an optional dynamic library name, for example
582: http://www.mcs.anl.gov/petsc/lib.a:coloringcreate()
584: Level: beginner
586: .seealso: MatGetColoring()
587: E*/
588: typedef char* MatColoringType;
589: #define MATCOLORING_NATURAL "natural"
590: #define MATCOLORING_SL "sl"
591: #define MATCOLORING_LF "lf"
592: #define MATCOLORING_ID "id"
594: EXTERN int MatGetColoring(Mat,MatColoringType,ISColoring*);
595: EXTERN int MatColoringRegister(char*,char*,char*,int(*)(Mat,MatColoringType,ISColoring *));
596: #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
597: #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,0)
598: #else
599: #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,d)
600: #endif
601: EXTERN int MatColoringRegisterAll(char *);
602: extern PetscTruth MatColoringRegisterAllCalled;
603: EXTERN int MatColoringRegisterDestroy(void);
604: EXTERN int MatColoringPatch(Mat,int,int,int *,ISColoring*);
606: /*S
607: MatFDColoring - Object for computing a sparse Jacobian via finite differences
608: and coloring
610: Level: beginner
612: Concepts: coloring, sparse Jacobian, finite differences
614: .seealso: MatFDColoringCreate()
615: S*/
616: typedef struct _p_MatFDColoring *MatFDColoring;
618: EXTERN int MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
619: EXTERN int MatFDColoringDestroy(MatFDColoring);
620: EXTERN int MatFDColoringView(MatFDColoring,PetscViewer);
621: EXTERN int MatFDColoringSetFunction(MatFDColoring,int (*)(void),void*);
622: EXTERN int MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
623: EXTERN int MatFDColoringSetFrequency(MatFDColoring,int);
624: EXTERN int MatFDColoringGetFrequency(MatFDColoring,int*);
625: EXTERN int MatFDColoringSetFromOptions(MatFDColoring);
626: EXTERN int MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *);
627: EXTERN int MatFDColoringApplyTS(Mat,MatFDColoring,PetscReal,Vec,MatStructure*,void *);
628: EXTERN int MatFDColoringSetRecompute(MatFDColoring);
629: EXTERN int MatFDColoringSetF(MatFDColoring,Vec);
631: /*
632: These routines are for partitioning matrices: currently used only
633: for adjacency matrix, MatCreateMPIAdj().
634: */
636: /*S
637: MatPartitioning - Object for managing the partitioning of a matrix or graph
639: Level: beginner
641: Concepts: partitioning
643: .seealso: MatParitioningCreate(), MatPartitioningType
644: S*/
645: typedef struct _p_MatPartitioning *MatPartitioning;
647: /*E
648: MatPartitioningType - String with the name of a PETSc matrix partitioing or the creation function
649: with an optional dynamic library name, for example
650: http://www.mcs.anl.gov/petsc/lib.a:partitioningcreate()
652: Level: beginner
654: .seealso: MatPartitioingCreate(), MatPartitioning
655: E*/
656: typedef char* MatPartitioningType;
657: #define MAT_PARTITIONING_CURRENT "current"
658: #define MAT_PARTITIONING_PARMETIS "parmetis"
660: EXTERN int MatPartitioningCreate(MPI_Comm,MatPartitioning*);
661: EXTERN int MatPartitioningSetType(MatPartitioning,MatPartitioningType);
662: EXTERN int MatPartitioningSetAdjacency(MatPartitioning,Mat);
663: EXTERN int MatPartitioningSetVertexWeights(MatPartitioning,int*);
664: EXTERN int MatPartitioningApply(MatPartitioning,IS*);
665: EXTERN int MatPartitioningDestroy(MatPartitioning);
667: EXTERN int MatPartitioningRegister(char*,char*,char*,int(*)(MatPartitioning));
668: #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
669: #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,0)
670: #else
671: #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,d)
672: #endif
674: EXTERN int MatPartitioningRegisterAll(char *);
675: extern PetscTruth MatPartitioningRegisterAllCalled;
676: EXTERN int MatPartitioningRegisterDestroy(void);
678: EXTERN int MatPartitioningView(MatPartitioning,PetscViewer);
679: EXTERN int MatPartitioningSetFromOptions(MatPartitioning);
680: EXTERN int MatPartitioningGetType(MatPartitioning,MatPartitioningType*);
682: EXTERN int MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
684: /*
685: If you add entries here you must also add them to finclude/petscmat.h
686: */
687: typedef enum { MATOP_SET_VALUES=0,
688: MATOP_GET_ROW=1,
689: MATOP_RESTORE_ROW=2,
690: MATOP_MULT=3,
691: MATOP_MULT_ADD=4,
692: MATOP_MULT_TRANSPOSE=5,
693: MATOP_MULT_TRANSPOSE_ADD=6,
694: MATOP_SOLVE=7,
695: MATOP_SOLVE_ADD=8,
696: MATOP_SOLVE_TRANSPOSE=9,
697: MATOP_SOLVE_TRANSPOSE_ADD=10,
698: MATOP_LUFACTOR=11,
699: MATOP_CHOLESKYFACTOR=12,
700: MATOP_RELAX=13,
701: MATOP_TRANSPOSE=14,
702: MATOP_GETINFO=15,
703: MATOP_EQUAL=16,
704: MATOP_GET_DIAGONAL=17,
705: MATOP_DIAGONAL_SCALE=18,
706: MATOP_NORM=19,
707: MATOP_ASSEMBLY_BEGIN=20,
708: MATOP_ASSEMBLY_END=21,
709: MATOP_COMPRESS=22,
710: MATOP_SET_OPTION=23,
711: MATOP_ZERO_ENTRIES=24,
712: MATOP_ZERO_ROWS=25,
713: MATOP_LUFACTOR_SYMBOLIC=26,
714: MATOP_LUFACTOR_NUMERIC=27,
715: MATOP_CHOLESKY_FACTOR_SYMBOLIC=28,
716: MATOP_CHOLESKY_FACTOR_NUMERIC=29,
717: MATOP_SETUP_PREALLOCATION=30,
718: MATOP_ILUFACTOR_SYMBOLIC=31,
719: MATOP_ICCFACTOR_SYMBOLIC=32,
720: MATOP_GET_ARRAY=33,
721: MATOP_RESTORE_ARRAY=34,
722: MATOP_DUPLCIATE=35,
723: MATOP_FORWARD_SOLVE=36,
724: MATOP_BACKWARD_SOLVE=37,
725: MATOP_ILUFACTOR=38,
726: MATOP_ICCFACTOR=39,
727: MATOP_AXPY=40,
728: MATOP_GET_SUBMATRICES=41,
729: MATOP_INCREASE_OVERLAP=42,
730: MATOP_GET_VALUES=43,
731: MATOP_COPY=44,
732: MATOP_PRINT_HELP=45,
733: MATOP_SCALE=46,
734: MATOP_SHIFT=47,
735: MATOP_DIAGONAL_SHIFT=48,
736: MATOP_ILUDT_FACTOR=49,
737: MATOP_GET_BLOCK_SIZE=50,
738: MATOP_GET_ROW_IJ=51,
739: MATOP_RESTORE_ROW_IJ=52,
740: MATOP_GET_COLUMN_IJ=53,
741: MATOP_RESTORE_COLUMN_IJ=54,
742: MATOP_FDCOLORING_CREATE=55,
743: MATOP_COLORING_PATCH=56,
744: MATOP_SET_UNFACTORED=57,
745: MATOP_PERMUTE=58,
746: MATOP_SET_VALUES_BLOCKED=59,
747: MATOP_GET_SUBMATRIX=60,
748: MATOP_DESTROY=61,
749: MATOP_VIEW=62,
750: MATOP_GET_MAPS=63,
751: MATOP_USE_SCALED_FORM=64,
752: MATOP_SCALE_SYSTEM=65,
753: MATOP_UNSCALE_SYSTEM=66,
754: MATOP_SET_LOCAL_TO_GLOBAL_MAPPING=67,
755: MATOP_SET_VALUES_LOCAL=68,
756: MATOP_ZERO_ROWS_LOCAL=69,
757: MATOP_GET_ROW_MAX=70,
758: MATOP_CONVERT=71,
759: MATOP_SET_COLORING=72,
760: MATOP_SET_VALUES_ADIC=73,
761: MATOP_SET_VALUES_ADIFOR=74,
762: MATOP_FD_COLORING_APPLY=75,
763: MATOP_SET_FROM_OPTIONS=76,
764: MATOP_MULT_CONSTRAINED=77,
765: MATOP_MULT_TRANSPOSE_CONSTRAINED=78,
766: MATOP_ILU_FACTOR_SYMBOLIC_CONSTRAINED=79,
767: MATOP_PERMUTE_SPARSIFY=80,
768: MATOP_MULT_MULTIPLE=81,
769: MATOP_SOLVE_MULTIPLE=82
770: } MatOperation;
771: EXTERN int MatHasOperation(Mat,MatOperation,PetscTruth*);
772: EXTERN int MatShellSetOperation(Mat,MatOperation,void(*)(void));
773: EXTERN int MatShellGetOperation(Mat,MatOperation,void(**)(void));
774: EXTERN int MatShellSetContext(Mat,void*);
776: /*
777: Codes for matrices stored on disk. By default they are
778: stored in a universal format. By changing the format with
779: PetscViewerSetFormat(viewer,PETSC_VIEWER_BINARY_NATIVE); the matrices will
780: be stored in a way natural for the matrix, for example dense matrices
781: would be stored as dense. Matrices stored this way may only be
782: read into matrices of the same time.
783: */
784: #define MATRIX_BINARY_FORMAT_DENSE -1
786: /*
787: New matrix classes not yet distributed
788: */
789: /*
790: MatAIJIndices is a data structure for storing the nonzero location information
791: for sparse matrices. Several matrices with identical nonzero structure can share
792: the same MatAIJIndices.
793: */
794: typedef struct _p_MatAIJIndices* MatAIJIndices;
796: EXTERN int MatCreateAIJIndices(int,int,int*,int*,PetscTruth,MatAIJIndices*);
797: EXTERN int MatCreateAIJIndicesEmpty(int,int,int*,PetscTruth,MatAIJIndices*);
798: EXTERN int MatAttachAIJIndices(MatAIJIndices,MatAIJIndices*);
799: EXTERN int MatDestroyAIJIndices(MatAIJIndices);
800: EXTERN int MatCopyAIJIndices(MatAIJIndices,MatAIJIndices*);
801: EXTERN int MatValidateAIJIndices(int,MatAIJIndices);
802: EXTERN int MatShiftAIJIndices(MatAIJIndices);
803: EXTERN int MatShrinkAIJIndices(MatAIJIndices);
804: EXTERN int MatTransposeAIJIndices(MatAIJIndices,MatAIJIndices*);
806: EXTERN int MatCreateSeqCSN(MPI_Comm,int,int,int*,int,Mat*);
807: EXTERN int MatCreateSeqCSN_Single(MPI_Comm,int,int,int*,int,Mat*);
808: EXTERN int MatCreateSeqCSNWithPrecision(MPI_Comm,int,int,int*,int,PetscScalarPrecision,Mat*);
810: EXTERN int MatCreateSeqCSNIndices(MPI_Comm,MatAIJIndices,int,Mat *);
811: EXTERN int MatCreateSeqCSNIndices_Single(MPI_Comm,MatAIJIndices,int,Mat *);
812: EXTERN int MatCreateSeqCSNIndicesWithPrecision(MPI_Comm,MatAIJIndices,int,PetscScalarPrecision,Mat *);
814: EXTERN int MatMPIBAIJSetHashTableFactor(Mat,PetscReal);
815: EXTERN int MatSeqAIJGetInodeSizes(Mat,int *,int *[],int *);
816: EXTERN int MatMPIRowbsGetColor(Mat,ISColoring *);
818: /*S
819: MatNullSpace - Object that removes a null space from a vector, i.e.
820: orthogonalizes the vector to a subsapce
822: Level: beginner
824: Concepts: matrix; linear operator, null space
826: Users manual sections:
827: . Section 4.15 Solving Singular Systems
829: .seealso: MatNullSpaceCreate()
830: S*/
831: typedef struct _p_MatNullSpace* MatNullSpace;
833: EXTERN int MatNullSpaceCreate(MPI_Comm,int,int,Vec *,MatNullSpace*);
834: EXTERN int MatNullSpaceDestroy(MatNullSpace);
835: EXTERN int MatNullSpaceRemove(MatNullSpace,Vec,Vec*);
836: EXTERN int MatNullSpaceAttach(Mat,MatNullSpace);
837: EXTERN int MatNullSpaceTest(MatNullSpace,Mat);
839: EXTERN int MatReorderingSeqSBAIJ(Mat A,IS isp);
840: EXTERN int MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
841: EXTERN int MatSeqSBAIJSetColumnIndices(Mat,int *);
844: EXTERN int MatCreateMAIJ(Mat,int,Mat*);
845: EXTERN int MatMAIJRedimension(Mat,int,Mat*);
846: EXTERN int MatMAIJGetAIJ(Mat,Mat*);
848: EXTERN int MatMPIAdjSetValues(Mat,int*,int*,int*);
850: EXTERN int MatComputeExplicitOperator(Mat,Mat*);
852: EXTERN int MatESISetType(Mat,char*);
853: EXTERN int MatESISetFromOptions(Mat);
855: EXTERN int MatMPIBAIJDiagonalScaleLocalSetUp(Mat,Vec);
856: EXTERN int MatMPIBAIJDiagonalScaleLocal(Mat,Vec);
858: EXTERN int PetscViewerMathematicaPutMatrix(PetscViewer, int, int, PetscReal *);
859: EXTERN int PetscViewerMathematicaPutCSRMatrix(PetscViewer, int, int, int *, int *, PetscReal *);
861: #endif