Actual source code: petscmat.h

petsc-main 2021-04-20
Report Typos and Errors
  1: /*
  2:      Include file for the matrix component of PETSc
  3: */
  4: #ifndef PETSCMAT_H
  5: #define PETSCMAT_H
  6: #include <petscvec.h>

  8: /*S
  9:      Mat - Abstract PETSc matrix object used to manage all linear operators in PETSc, even those without
 10:            an explicit sparse representation (such as matrix-free operators)

 12:    Level: beginner

 14: .seealso:  MatCreate(), MatType, MatSetType(), MatDestroy()
 15: S*/
 16: typedef struct _p_Mat*           Mat;

 18: /*J
 19:     MatType - String with the name of a PETSc matrix type

 21:    Level: beginner

 23: .seealso: MatSetType(), Mat, MatSolverType, MatRegister()
 24: J*/
 25: typedef const char* MatType;
 26: #define MATSAME            "same"
 27: #define MATMAIJ            "maij"
 28: #define MATSEQMAIJ         "seqmaij"
 29: #define MATMPIMAIJ         "mpimaij"
 30: #define MATKAIJ            "kaij"
 31: #define MATSEQKAIJ         "seqkaij"
 32: #define MATMPIKAIJ         "mpikaij"
 33: #define MATIS              "is"
 34: #define MATAIJ             "aij"
 35: #define MATSEQAIJ          "seqaij"
 36: #define MATMPIAIJ          "mpiaij"
 37: #define MATAIJCRL          "aijcrl"
 38: #define MATSEQAIJCRL       "seqaijcrl"
 39: #define MATMPIAIJCRL       "mpiaijcrl"
 40: #define MATAIJCUSPARSE     "aijcusparse"
 41: #define MATSEQAIJCUSPARSE  "seqaijcusparse"
 42: #define MATMPIAIJCUSPARSE  "mpiaijcusparse"
 43: #define MATAIJKOKKOS       "aijkokkos"
 44: #define MATSEQAIJKOKKOS    "seqaijkokkos"
 45: #define MATMPIAIJKOKKOS    "mpiaijkokkos"
 46: #define MATAIJVIENNACL     "aijviennacl"
 47: #define MATSEQAIJVIENNACL  "seqaijviennacl"
 48: #define MATMPIAIJVIENNACL  "mpiaijviennacl"
 49: #define MATAIJPERM         "aijperm"
 50: #define MATSEQAIJPERM      "seqaijperm"
 51: #define MATMPIAIJPERM      "mpiaijperm"
 52: #define MATAIJSELL         "aijsell"
 53: #define MATSEQAIJSELL      "seqaijsell"
 54: #define MATMPIAIJSELL      "mpiaijsell"
 55: #define MATAIJMKL          "aijmkl"
 56: #define MATSEQAIJMKL       "seqaijmkl"
 57: #define MATMPIAIJMKL       "mpiaijmkl"
 58: #define MATBAIJMKL         "baijmkl"
 59: #define MATSEQBAIJMKL      "seqbaijmkl"
 60: #define MATMPIBAIJMKL      "mpibaijmkl"
 61: #define MATSHELL           "shell"
 62: #define MATDENSE           "dense"
 63: #define MATDENSECUDA       "densecuda"
 64: #define MATSEQDENSE        "seqdense"
 65: #define MATSEQDENSECUDA    "seqdensecuda"
 66: #define MATMPIDENSE        "mpidense"
 67: #define MATMPIDENSECUDA    "mpidensecuda"
 68: #define MATELEMENTAL       "elemental"
 69: #define MATSCALAPACK       "scalapack"
 70: #define MATBAIJ            "baij"
 71: #define MATSEQBAIJ         "seqbaij"
 72: #define MATMPIBAIJ         "mpibaij"
 73: #define MATMPIADJ          "mpiadj"
 74: #define MATSBAIJ           "sbaij"
 75: #define MATSEQSBAIJ        "seqsbaij"
 76: #define MATMPISBAIJ        "mpisbaij"
 77: #define MATMFFD            "mffd"
 78: #define MATNORMAL          "normal"
 79: #define MATNORMALHERMITIAN "normalh"
 80: #define MATLRC             "lrc"
 81: #define MATSCATTER         "scatter"
 82: #define MATBLOCKMAT        "blockmat"
 83: #define MATCOMPOSITE       "composite"
 84: #define MATFFT             "fft"
 85: #define MATFFTW            "fftw"
 86: #define MATSEQCUFFT        "seqcufft"
 87: #define MATTRANSPOSEMAT    "transpose"
 88: #define MATSCHURCOMPLEMENT "schurcomplement"
 89: #define MATPYTHON          "python"
 90: #define MATHYPRE           "hypre"
 91: #define MATHYPRESTRUCT     "hyprestruct"
 92: #define MATHYPRESSTRUCT    "hypresstruct"
 93: #define MATSUBMATRIX       "submatrix"
 94: #define MATLOCALREF        "localref"
 95: #define MATNEST            "nest"
 96: #define MATPREALLOCATOR    "preallocator"
 97: #define MATSELL            "sell"
 98: #define MATSEQSELL         "seqsell"
 99: #define MATMPISELL         "mpisell"
100: #define MATDUMMY           "dummy"
101: #define MATLMVM            "lmvm"
102: #define MATLMVMDFP         "lmvmdfp"
103: #define MATLMVMBFGS        "lmvmbfgs"
104: #define MATLMVMSR1         "lmvmsr1"
105: #define MATLMVMBROYDEN     "lmvmbroyden"
106: #define MATLMVMBADBROYDEN  "lmvmbadbroyden"
107: #define MATLMVMSYMBROYDEN  "lmvmsymbroyden"
108: #define MATLMVMSYMBADBROYDEN "lmvmsymbadbroyden"
109: #define MATLMVMDIAGBROYDEN   "lmvmdiagbroyden"
110: #define MATCONSTANTDIAGONAL  "constantdiagonal"
111: #define MATHARA              "hara"

113: /*J
114:     MatSolverType - String with the name of a PETSc matrix solver type.

116:     For example: "petsc" indicates what PETSc provides, "superlu_dist" the parallel SuperLU_DIST package etc

118:    Level: beginner

120:    Notes:  MATSOLVERUMFPACK, MATSOLVERCHOLMOD, MATSOLVERKLU, MATSOLVERSPQR form the SuiteSparse package for which you can use --download-suitesparse

122: .seealso: MatGetFactor(), PCFactorSetMatSolverType(), PCFactorGetMatSolverType()
123: J*/
124: typedef const char* MatSolverType;
125: #define MATSOLVERSUPERLU          "superlu"
126: #define MATSOLVERSUPERLU_DIST     "superlu_dist"
127: #define MATSOLVERSTRUMPACK        "strumpack"
128: #define MATSOLVERUMFPACK          "umfpack"
129: #define MATSOLVERCHOLMOD          "cholmod"
130: #define MATSOLVERKLU              "klu"
131: #define MATSOLVERSPARSEELEMENTAL  "sparseelemental"
132: #define MATSOLVERELEMENTAL        "elemental"
133: #define MATSOLVERSCALAPACK        "scalapack"
134: #define MATSOLVERESSL             "essl"
135: #define MATSOLVERLUSOL            "lusol"
136: #define MATSOLVERMUMPS            "mumps"
137: #define MATSOLVERMKL_PARDISO      "mkl_pardiso"
138: #define MATSOLVERMKL_CPARDISO     "mkl_cpardiso"
139: #define MATSOLVERPASTIX           "pastix"
140: #define MATSOLVERMATLAB           "matlab"
141: #define MATSOLVERPETSC            "petsc"
142: #define MATSOLVERBAS              "bas"
143: #define MATSOLVERCUSPARSE         "cusparse"
144: #define MATSOLVERCUSPARSEBAND     "cusparseband"
145: #define MATSOLVERCUDA             "cuda"
146: #define MATSOLVERKOKKOS           "kokkos"
147: #define MATSOLVERKOKKOSDEVICE     "kokkosdevice"

149: /*E
150:     MatFactorType - indicates what type of factorization is requested

152:     Level: beginner

154:    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h

156: .seealso: MatSolverType, MatGetFactor(), MatGetFactorAvailable(), MatSolverTypeRegister()
157: E*/
158: typedef enum {MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC, MAT_FACTOR_ILUDT, MAT_FACTOR_QR, MAT_FACTOR_NUM_TYPES} MatFactorType;
159: PETSC_EXTERN const char *const MatFactorTypes[];

161: PETSC_EXTERN PetscErrorCode MatGetFactor(Mat,MatSolverType,MatFactorType,Mat*);
162: PETSC_EXTERN PetscErrorCode MatGetFactorAvailable(Mat,MatSolverType,MatFactorType,PetscBool *);
163: PETSC_EXTERN PetscErrorCode MatFactorGetCanUseOrdering(Mat, PetscBool*);
164: PETSC_DEPRECATED_FUNCTION("Use MatFactorGetCanUseOrdering() (since version 3.15)") PETSC_STATIC_INLINE PetscErrorCode MatFactorGetUseOrdering(Mat A,PetscBool *b) {return MatFactorGetCanUseOrdering(A,b);}
165: PETSC_EXTERN PetscErrorCode MatFactorGetSolverType(Mat,MatSolverType*);
166: PETSC_EXTERN PetscErrorCode MatGetFactorType(Mat,MatFactorType*);
167: PETSC_EXTERN PetscErrorCode MatSetFactorType(Mat,MatFactorType);
168: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*MatSolverFunction)(Mat,MatFactorType,Mat*);
169: PETSC_EXTERN PetscErrorCode MatSolverTypeRegister(MatSolverType,MatType,MatFactorType,MatSolverFunction);
170: PETSC_EXTERN PetscErrorCode MatSolverTypeGet(MatSolverType,MatType,MatFactorType,PetscBool*,PetscBool*,MatSolverFunction*);
171: typedef MatSolverType MatSolverPackage PETSC_DEPRECATED_TYPEDEF("Use MatSolverType (since version 3.9)");
172: PETSC_DEPRECATED_FUNCTION("Use MatSolverTypeRegister() (since version 3.9)") PETSC_STATIC_INLINE PetscErrorCode MatSolverPackageRegister(MatSolverType stype,MatType mtype,MatFactorType ftype,MatSolverFunction f)
173: { return MatSolverTypeRegister(stype,mtype,ftype,f); }
174: PETSC_DEPRECATED_FUNCTION("Use MatSolverTypeGet() (since version 3.9)") PETSC_STATIC_INLINE PetscErrorCode MatSolverPackageGet(MatSolverType stype,MatType mtype,MatFactorType ftype,PetscBool *foundmtype,PetscBool *foundstype,MatSolverFunction *f)
175: { return MatSolverTypeGet(stype,mtype,ftype,foundmtype,foundstype,f); }

177: /*E
178:     MatProductType - indicates what type of matrix product is requested

180:     Level: beginner

182: .seealso: MatSolverType, MatProductSetType()
183: E*/
184: typedef enum {MATPRODUCT_UNSPECIFIED=0,MATPRODUCT_AB,MATPRODUCT_AtB,MATPRODUCT_ABt,MATPRODUCT_PtAP,MATPRODUCT_RARt,MATPRODUCT_ABC} MatProductType;
185: PETSC_EXTERN const char *const MatProductTypes[];

187: /*J
188:     MatProductAlgorithm - String with the name of an algorithm for a PETSc matrix product implementation

190:    Level: beginner

192: .seealso: MatSetType(), Mat, MatSolverType, MatRegister(), MatProductSetAlgorithm(), MatProductType
193: J*/
194: typedef const char* MatProductAlgorithm;
195: #define MATPRODUCTALGORITHM_DEFAULT "default"

197: PETSC_EXTERN PetscErrorCode MatProductCreate(Mat,Mat,Mat,Mat*);
198: PETSC_EXTERN PetscErrorCode MatProductCreateWithMat(Mat,Mat,Mat,Mat);
199: PETSC_EXTERN PetscErrorCode MatProductSetType(Mat,MatProductType);
200: PETSC_EXTERN PetscErrorCode MatProductSetAlgorithm(Mat,MatProductAlgorithm);
201: PETSC_EXTERN PetscErrorCode MatProductSetFill(Mat,PetscReal);
202: PETSC_EXTERN PetscErrorCode MatProductSetFromOptions(Mat);
203: PETSC_EXTERN PetscErrorCode MatProductSymbolic(Mat);
204: PETSC_EXTERN PetscErrorCode MatProductNumeric(Mat);
205: PETSC_EXTERN PetscErrorCode MatProductReplaceMats(Mat,Mat,Mat,Mat);
206: PETSC_EXTERN PetscErrorCode MatProductClear(Mat);
207: PETSC_EXTERN PetscErrorCode MatProductView(Mat,PetscViewer);

209: /* Logging support */
210: #define    MAT_FILE_CLASSID 1211216    /* used to indicate matrices in binary files */
211: PETSC_EXTERN PetscClassId MAT_CLASSID;
212: PETSC_EXTERN PetscClassId MAT_COLORING_CLASSID;
213: PETSC_EXTERN PetscClassId MAT_FDCOLORING_CLASSID;
214: PETSC_EXTERN PetscClassId MAT_TRANSPOSECOLORING_CLASSID;
215: PETSC_EXTERN PetscClassId MAT_PARTITIONING_CLASSID;
216: PETSC_EXTERN PetscClassId MAT_COARSEN_CLASSID;
217: PETSC_EXTERN PetscClassId MAT_NULLSPACE_CLASSID;
218: PETSC_EXTERN PetscClassId MATMFFD_CLASSID;

220: /*E
221:     MatReuse - Indicates if matrices obtained from a previous call to MatCreateSubMatrices(), MatCreateSubMatrix(), MatConvert() or several other functions
222:      are to be reused to store the new matrix values.

224: $  MAT_INITIAL_MATRIX - create a new matrix
225: $  MAT_REUSE_MATRIX - reuse the matrix created with a previous call that used MAT_INITIAL_MATRIX
226: $  MAT_INPLACE_MATRIX - replace the first input matrix with the new matrix (not applicable to all functions)
227: $  MAT_IGNORE_MATRIX - do not create a new matrix or reuse a give matrix, just ignore that matrix argument (not applicable to all functions)

229:     Level: beginner

231:    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h

233: .seealso: MatCreateSubMatrices(), MatCreateSubMatrix(), MatDestroyMatrices(), MatConvert()
234: E*/
235: typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX,MAT_IGNORE_MATRIX,MAT_INPLACE_MATRIX} MatReuse;

237: /*E
238:     MatCreateSubMatrixOption - Indicates if matrices obtained from a call to MatCreateSubMatrices()
239:      include the matrix values. Currently it is only used by MatGetSeqNonzeroStructure().

241:     Level: beginner

243: .seealso: MatGetSeqNonzeroStructure()
244: E*/
245: typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatCreateSubMatrixOption;

247: PETSC_EXTERN PetscErrorCode MatInitializePackage(void);

249: PETSC_EXTERN PetscErrorCode MatCreate(MPI_Comm,Mat*);
250: PETSC_EXTERN PetscErrorCode MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt);
251: PETSC_EXTERN PetscErrorCode MatSetType(Mat,MatType);
252: PETSC_EXTERN PetscErrorCode MatGetVecType(Mat,VecType*);
253: PETSC_EXTERN PetscErrorCode MatSetVecType(Mat,VecType);
254: PETSC_EXTERN PetscErrorCode MatSetFromOptions(Mat);
255: PETSC_EXTERN PetscErrorCode MatViewFromOptions(Mat,PetscObject,const char[]);
256: PETSC_EXTERN PetscErrorCode MatRegister(const char[],PetscErrorCode(*)(Mat));
257: PETSC_EXTERN PetscErrorCode MatRegisterRootName(const char[],const char[],const char[]);
258: PETSC_EXTERN PetscErrorCode MatSetOptionsPrefix(Mat,const char[]);
259: PETSC_EXTERN PetscErrorCode MatAppendOptionsPrefix(Mat,const char[]);
260: PETSC_EXTERN PetscErrorCode MatGetOptionsPrefix(Mat,const char*[]);
261: PETSC_EXTERN PetscErrorCode MatSetErrorIfFailure(Mat,PetscBool);

263: PETSC_EXTERN PetscFunctionList MatList;
264: PETSC_EXTERN PetscFunctionList MatColoringList;
265: PETSC_EXTERN PetscFunctionList MatPartitioningList;

267: /*E
268:     MatStructure - Indicates if two matrices have the same nonzero structure

270:     Level: beginner

272: $  SAME_NONZERO_PATTERN  - the two matrices have identical nonzero patterns
273: $  DIFFERENT_NONZERO_PATTERN - the two matrices may have different nonzero patterns
274: $  SUBSET_NONZERO_PATTERN - the nonzero pattern of the second matrix is a subset of the nonzero pattern of the first matrix
275: $  UNKNOWN_NONZERO_PATTERN - there is no known relationship between the nonzero patterns. In this case the implementations may try to detect a relationship to optimize the operation

277:    Developer Notes:
278:      Any additions/changes here MUST also be made in src/mat/f90-mod/petscmat.h

280: .seealso: MatCopy(), MatAXPY(), MatAYPX()
281: E*/
282: typedef enum {DIFFERENT_NONZERO_PATTERN,SUBSET_NONZERO_PATTERN,SAME_NONZERO_PATTERN,UNKNOWN_NONZERO_PATTERN} MatStructure;
283: PETSC_EXTERN const char *const MatStructures[];

285: #if defined PETSC_HAVE_MKL_SPARSE
286: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJMKL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
287: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJMKL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
288: PETSC_EXTERN PetscErrorCode MatCreateBAIJMKL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
289: PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJMKL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
290: #endif

292: PETSC_EXTERN PetscErrorCode MatCreateSeqSELL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
293: PETSC_EXTERN PetscErrorCode MatCreateSELL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
294: PETSC_EXTERN PetscErrorCode MatSeqSELLSetPreallocation(Mat,PetscInt,const PetscInt[]);
295: PETSC_EXTERN PetscErrorCode MatMPISELLSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);

297: PETSC_EXTERN PetscErrorCode MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
298: PETSC_EXTERN PetscErrorCode MatCreateDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
299: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
300: PETSC_EXTERN PetscErrorCode MatCreateAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
301: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat*);
302: PETSC_EXTERN PetscErrorCode MatUpdateMPIAIJWithArrays(Mat,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
303: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*);
304: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithSeqAIJ(MPI_Comm,Mat,Mat,const PetscInt[],Mat*);

306: PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
307: PETSC_EXTERN PetscErrorCode MatCreateBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
308: PETSC_EXTERN PetscErrorCode MatCreateMPIBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat*);

310: PETSC_EXTERN PetscErrorCode MatSetPreallocationCOO(Mat,PetscInt,const PetscInt[],const PetscInt[]);
311: PETSC_EXTERN PetscErrorCode MatSetValuesCOO(Mat,const PetscScalar[],InsertMode);

313: PETSC_EXTERN PetscErrorCode MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*);
314: PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);

316: PETSC_EXTERN PetscErrorCode MatCreateSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
317: PETSC_EXTERN PetscErrorCode MatCreateMPISBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
318: PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
319: PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
320: PETSC_EXTERN PetscErrorCode MatXAIJSetPreallocation(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscInt[],const PetscInt[]);

322: PETSC_EXTERN PetscErrorCode MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*);
323: PETSC_EXTERN PetscErrorCode MatCreateNormal(Mat,Mat*);
324: PETSC_EXTERN PetscErrorCode MatCreateNormalHermitian(Mat,Mat*);
325: PETSC_EXTERN PetscErrorCode MatCreateLRC(Mat,Mat,Vec,Mat,Mat*);
326: PETSC_EXTERN PetscErrorCode MatLRCGetMats(Mat,Mat*,Mat*,Vec*,Mat*);
327: PETSC_EXTERN PetscErrorCode MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,ISLocalToGlobalMapping,Mat*);
328: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
329: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);

331: PETSC_EXTERN PetscErrorCode MatCreateScatter(MPI_Comm,VecScatter,Mat*);
332: PETSC_EXTERN PetscErrorCode MatScatterSetVecScatter(Mat,VecScatter);
333: PETSC_EXTERN PetscErrorCode MatScatterGetVecScatter(Mat,VecScatter*);
334: PETSC_EXTERN PetscErrorCode MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*);
335: PETSC_EXTERN PetscErrorCode MatCompositeAddMat(Mat,Mat);
336: PETSC_EXTERN PetscErrorCode MatCompositeMerge(Mat);
337: typedef enum {MAT_COMPOSITE_MERGE_RIGHT,MAT_COMPOSITE_MERGE_LEFT} MatCompositeMergeType;
338: PETSC_EXTERN PetscErrorCode MatCompositeSetMergeType(Mat,MatCompositeMergeType);
339: PETSC_EXTERN PetscErrorCode MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*);
340: typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType;
341: PETSC_EXTERN PetscErrorCode MatCompositeSetType(Mat,MatCompositeType);
342: PETSC_EXTERN PetscErrorCode MatCompositeGetType(Mat,MatCompositeType*);
343: PETSC_EXTERN PetscErrorCode MatCompositeSetMatStructure(Mat,MatStructure);
344: PETSC_EXTERN PetscErrorCode MatCompositeGetMatStructure(Mat,MatStructure*);
345: PETSC_EXTERN PetscErrorCode MatCompositeGetNumberMat(Mat,PetscInt*);
346: PETSC_EXTERN PetscErrorCode MatCompositeGetMat(Mat,PetscInt,Mat*);
347: PETSC_EXTERN PetscErrorCode MatCompositeSetScalings(Mat,const PetscScalar*);

349: PETSC_EXTERN PetscErrorCode MatCreateFFT(MPI_Comm,PetscInt,const PetscInt[],MatType,Mat*);
350: PETSC_EXTERN PetscErrorCode MatCreateSeqCUFFT(MPI_Comm,PetscInt,const PetscInt[],Mat*);

352: PETSC_EXTERN PetscErrorCode MatCreateTranspose(Mat,Mat*);
353: PETSC_EXTERN PetscErrorCode MatTransposeGetMat(Mat,Mat*);
354: PETSC_EXTERN PetscErrorCode MatCreateHermitianTranspose(Mat,Mat*);
355: PETSC_EXTERN PetscErrorCode MatHermitianTransposeGetMat(Mat,Mat*);
356: PETSC_EXTERN PetscErrorCode MatCreateSubMatrixVirtual(Mat,IS,IS,Mat*);
357: PETSC_EXTERN PetscErrorCode MatSubMatrixVirtualUpdate(Mat,Mat,IS,IS);
358: PETSC_EXTERN PetscErrorCode MatCreateLocalRef(Mat,IS,IS,Mat*);
359: PETSC_EXTERN PetscErrorCode MatCreateConstantDiagonal(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar,Mat*);

361: #if defined(PETSC_HAVE_HYPRE)
362: PETSC_EXTERN PetscErrorCode MatHYPRESetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
363: #endif

365: PETSC_EXTERN PetscErrorCode MatPythonSetType(Mat,const char[]);

367: PETSC_EXTERN PetscErrorCode MatResetPreallocation(Mat);
368: PETSC_EXTERN PetscErrorCode MatSetUp(Mat);
369: PETSC_EXTERN PetscErrorCode MatDestroy(Mat*);
370: PETSC_EXTERN PetscErrorCode MatGetNonzeroState(Mat,PetscObjectState*);

372: PETSC_EXTERN PetscErrorCode MatConjugate(Mat);
373: PETSC_EXTERN PetscErrorCode MatRealPart(Mat);
374: PETSC_EXTERN PetscErrorCode MatImaginaryPart(Mat);
375: PETSC_EXTERN PetscErrorCode MatGetDiagonalBlock(Mat,Mat*);
376: PETSC_EXTERN PetscErrorCode MatGetTrace(Mat,PetscScalar*);
377: PETSC_EXTERN PetscErrorCode MatInvertBlockDiagonal(Mat,const PetscScalar **);
378: PETSC_EXTERN PetscErrorCode MatInvertVariableBlockDiagonal(Mat,PetscInt,const PetscInt*,PetscScalar*);
379: PETSC_EXTERN PetscErrorCode MatInvertBlockDiagonalMat(Mat,Mat);

381: /* ------------------------------------------------------------*/
382: PETSC_EXTERN PetscErrorCode MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
383: PETSC_EXTERN PetscErrorCode MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
384: PETSC_EXTERN PetscErrorCode MatSetValuesRow(Mat,PetscInt,const PetscScalar[]);
385: PETSC_EXTERN PetscErrorCode MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]);
386: PETSC_EXTERN PetscErrorCode MatSetValuesBatch(Mat,PetscInt,PetscInt,PetscInt[],const PetscScalar[]);
387: PETSC_EXTERN PetscErrorCode MatSetRandom(Mat,PetscRandom);

389: /*S
390:      MatStencil - Data structure (C struct) for storing information about a single row or
391:         column of a matrix as indexed on an associated grid. These are arguments to MatSetStencil() and MatSetBlockStencil()

393:    The i,j, and k represent the logical coordinates over the entire grid (for 2 and 1 dimensional problems the k and j entries are ignored).
394:    The c represents the the degrees of freedom at each grid point (the dof argument to DMDASetDOF()). If dof is 1 then this entry is ignored.

396:    For stencil access to vectors see DMDAVecGetArray(), DMDAVecGetArrayF90().

398:    Fortran usage is different, see MatSetValuesStencil() for details.

400:    Level: beginner

402: .seealso:  MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockedStencil(), DMDAVecGetArray(), DMDAVecGetArrayF90()
403: S*/
404: typedef struct {
405:   PetscInt k,j,i,c;
406: } MatStencil;

408: PETSC_EXTERN PetscErrorCode MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
409: PETSC_EXTERN PetscErrorCode MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
410: PETSC_EXTERN PetscErrorCode MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt);

412: /*E
413:     MatAssemblyType - Indicates if the matrix is now to be used, or if you plan
414:      to continue to add or insert values to it

416:     Level: beginner

418: .seealso: MatAssemblyBegin(), MatAssemblyEnd()
419: E*/
420: typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType;
421: PETSC_EXTERN PetscErrorCode MatAssemblyBegin(Mat,MatAssemblyType);
422: PETSC_EXTERN PetscErrorCode MatAssemblyEnd(Mat,MatAssemblyType);
423: PETSC_EXTERN PetscErrorCode MatAssembled(Mat,PetscBool *);



427: /*E
428:     MatOption - Options that may be set for a matrix and its behavior or storage

430:     Level: beginner

432:    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
433:    Any additions/changes here must also be made in src/mat/interface/dlregismat.c in MatOptions[]

435:    Developer Notes:
436:     Entries that are negative need not be called collectively by all processes.

438: .seealso: MatSetOption()
439: E*/
440: typedef enum {MAT_OPTION_MIN = -3,
441:               MAT_UNUSED_NONZERO_LOCATION_ERR = -2,
442:               MAT_ROW_ORIENTED = -1,
443:               MAT_SYMMETRIC = 1,
444:               MAT_STRUCTURALLY_SYMMETRIC = 2,
445:               MAT_FORCE_DIAGONAL_ENTRIES = 3,
446:               MAT_IGNORE_OFF_PROC_ENTRIES = 4,
447:               MAT_USE_HASH_TABLE = 5,
448:               MAT_KEEP_NONZERO_PATTERN = 6,
449:               MAT_IGNORE_ZERO_ENTRIES = 7,
450:               MAT_USE_INODES = 8,
451:               MAT_HERMITIAN = 9,
452:               MAT_SYMMETRY_ETERNAL = 10,
453:               MAT_NEW_NONZERO_LOCATION_ERR = 11,
454:               MAT_IGNORE_LOWER_TRIANGULAR = 12,
455:               MAT_ERROR_LOWER_TRIANGULAR = 13,
456:               MAT_GETROW_UPPERTRIANGULAR = 14,
457:               MAT_SPD = 15,
458:               MAT_NO_OFF_PROC_ZERO_ROWS = 16,
459:               MAT_NO_OFF_PROC_ENTRIES = 17,
460:               MAT_NEW_NONZERO_LOCATIONS = 18,
461:               MAT_NEW_NONZERO_ALLOCATION_ERR = 19,
462:               MAT_SUBSET_OFF_PROC_ENTRIES = 20,
463:               MAT_SUBMAT_SINGLEIS = 21,
464:               MAT_STRUCTURE_ONLY = 22,
465:               MAT_SORTED_FULL = 23,
466:               MAT_FORM_EXPLICIT_TRANSPOSE = 24,
467:               MAT_OPTION_MAX = 25} MatOption;

469: PETSC_EXTERN const char *const *MatOptions;
470: PETSC_EXTERN PetscErrorCode MatSetOption(Mat,MatOption,PetscBool);
471: PETSC_EXTERN PetscErrorCode MatGetOption(Mat,MatOption,PetscBool*);
472: PETSC_EXTERN PetscErrorCode MatPropagateSymmetryOptions(Mat,Mat);
473: PETSC_EXTERN PetscErrorCode MatGetType(Mat,MatType*);

475: PETSC_EXTERN PetscErrorCode MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
476: PETSC_EXTERN PetscErrorCode MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
477: PETSC_EXTERN PetscErrorCode MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
478: PETSC_EXTERN PetscErrorCode MatGetRowUpperTriangular(Mat);
479: PETSC_EXTERN PetscErrorCode MatRestoreRowUpperTriangular(Mat);
480: PETSC_EXTERN PetscErrorCode MatGetColumnVector(Mat,Vec,PetscInt);
481: PETSC_EXTERN PetscErrorCode MatSeqAIJGetArray(Mat,PetscScalar *[]);
482: PETSC_EXTERN PetscErrorCode MatSeqAIJGetArrayRead(Mat,const PetscScalar *[]);
483: PETSC_EXTERN PetscErrorCode MatSeqAIJRestoreArray(Mat,PetscScalar *[]);
484: PETSC_EXTERN PetscErrorCode MatSeqAIJRestoreArrayRead(Mat,const PetscScalar *[]);
485: PETSC_EXTERN PetscErrorCode MatSeqAIJGetMaxRowNonzeros(Mat,PetscInt*);
486: PETSC_EXTERN PetscErrorCode MatSeqAIJSetValuesLocalFast(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
487: PETSC_EXTERN PetscErrorCode MatSeqAIJSetType(Mat,MatType);
488: PETSC_EXTERN PetscErrorCode MatSeqAIJRegister(const char[],PetscErrorCode (*)(Mat,MatType,MatReuse,Mat *));
489: PETSC_EXTERN PetscFunctionList MatSeqAIJList;
490: PETSC_EXTERN PetscErrorCode MatSeqBAIJGetArray(Mat,PetscScalar *[]);
491: PETSC_EXTERN PetscErrorCode MatSeqBAIJRestoreArray(Mat,PetscScalar *[]);
492: PETSC_EXTERN PetscErrorCode MatSeqSBAIJGetArray(Mat,PetscScalar *[]);
493: PETSC_EXTERN PetscErrorCode MatSeqSBAIJRestoreArray(Mat,PetscScalar *[]);
494: PETSC_EXTERN PetscErrorCode MatDenseGetArray(Mat,PetscScalar *[]);
495: PETSC_EXTERN PetscErrorCode MatDenseRestoreArray(Mat,PetscScalar *[]);
496: PETSC_EXTERN PetscErrorCode MatDensePlaceArray(Mat,const PetscScalar[]);
497: PETSC_EXTERN PetscErrorCode MatDenseReplaceArray(Mat,const PetscScalar[]);
498: PETSC_EXTERN PetscErrorCode MatDenseResetArray(Mat);
499: PETSC_EXTERN PetscErrorCode MatDenseGetArrayRead(Mat,const PetscScalar *[]);
500: PETSC_EXTERN PetscErrorCode MatDenseRestoreArrayRead(Mat,const PetscScalar *[]);
501: PETSC_EXTERN PetscErrorCode MatDenseGetArrayWrite(Mat,PetscScalar *[]);
502: PETSC_EXTERN PetscErrorCode MatDenseRestoreArrayWrite(Mat,PetscScalar *[]);
503: PETSC_EXTERN PetscErrorCode MatGetBlockSize(Mat,PetscInt *);
504: PETSC_EXTERN PetscErrorCode MatSetBlockSize(Mat,PetscInt);
505: PETSC_EXTERN PetscErrorCode MatGetBlockSizes(Mat,PetscInt *,PetscInt *);
506: PETSC_EXTERN PetscErrorCode MatSetBlockSizes(Mat,PetscInt,PetscInt);
507: PETSC_EXTERN PetscErrorCode MatSetBlockSizesFromMats(Mat,Mat,Mat);
508: PETSC_EXTERN PetscErrorCode MatSetVariableBlockSizes(Mat,PetscInt,PetscInt*);
509: PETSC_EXTERN PetscErrorCode MatGetVariableBlockSizes(Mat,PetscInt*,const PetscInt**);

511: PETSC_EXTERN PetscErrorCode MatDenseGetColumn(Mat,PetscInt,PetscScalar*[]);
512: PETSC_EXTERN PetscErrorCode MatDenseRestoreColumn(Mat,PetscScalar*[]);
513: PETSC_EXTERN PetscErrorCode MatDenseGetColumnVec(Mat,PetscInt,Vec*);
514: PETSC_EXTERN PetscErrorCode MatDenseRestoreColumnVec(Mat,PetscInt,Vec*);
515: PETSC_EXTERN PetscErrorCode MatDenseGetColumnVecRead(Mat,PetscInt,Vec*);
516: PETSC_EXTERN PetscErrorCode MatDenseRestoreColumnVecRead(Mat,PetscInt,Vec*);
517: PETSC_EXTERN PetscErrorCode MatDenseGetColumnVecWrite(Mat,PetscInt,Vec*);
518: PETSC_EXTERN PetscErrorCode MatDenseRestoreColumnVecWrite(Mat,PetscInt,Vec*);
519: PETSC_EXTERN PetscErrorCode MatDenseGetSubMatrix(Mat,PetscInt,PetscInt,Mat*);
520: PETSC_EXTERN PetscErrorCode MatDenseRestoreSubMatrix(Mat,Mat*);

522: PETSC_EXTERN PetscErrorCode MatMult(Mat,Vec,Vec);
523: PETSC_EXTERN PetscErrorCode MatMultDiagonalBlock(Mat,Vec,Vec);
524: PETSC_EXTERN PetscErrorCode MatMultAdd(Mat,Vec,Vec,Vec);
525: PETSC_EXTERN PetscErrorCode MatMultTranspose(Mat,Vec,Vec);
526: PETSC_EXTERN PetscErrorCode MatMultHermitianTranspose(Mat,Vec,Vec);
527: PETSC_EXTERN PetscErrorCode MatIsTranspose(Mat,Mat,PetscReal,PetscBool *);
528: PETSC_EXTERN PetscErrorCode MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *);
529: PETSC_EXTERN PetscErrorCode MatMultTransposeAdd(Mat,Vec,Vec,Vec);
530: PETSC_EXTERN PetscErrorCode MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec);
531: PETSC_EXTERN PetscErrorCode MatMultConstrained(Mat,Vec,Vec);
532: PETSC_EXTERN PetscErrorCode MatMultTransposeConstrained(Mat,Vec,Vec);
533: PETSC_EXTERN PetscErrorCode MatMatSolve(Mat,Mat,Mat);
534: PETSC_EXTERN PetscErrorCode MatMatSolveTranspose(Mat,Mat,Mat);
535: PETSC_EXTERN PetscErrorCode MatMatTransposeSolve(Mat,Mat,Mat);
536: PETSC_EXTERN PetscErrorCode MatResidual(Mat,Vec,Vec,Vec);

538: /*E
539:     MatDuplicateOption - Indicates if a duplicated sparse matrix should have
540:   its numerical values copied over or just its nonzero structure.

542:     Level: beginner

544:    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h

546: $   MAT_DO_NOT_COPY_VALUES    - Create a matrix using the same nonzero pattern as the original matrix,
547: $                               with zeros for the numerical values.
548: $   MAT_COPY_VALUES           - Create a matrix with the same nonzero pattern as the original matrix
549: $                               and with the same numerical values.
550: $   MAT_SHARE_NONZERO_PATTERN - Create a matrix that shares the nonzero structure with the previous matrix
551: $                               and does not copy it, using zeros for the numerical values. The parent and
552: $                               child matrices will share their index (i and j) arrays, and you cannot
553: $                               insert new nonzero entries into either matrix.

555: Notes:
556:     Many matrix types (including SeqAIJ) do not support the MAT_SHARE_NONZERO_PATTERN optimization; in
557: this case the behavior is as if MAT_DO_NOT_COPY_VALUES has been specified.

559: .seealso: MatDuplicate()
560: E*/
561: typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption;

563: PETSC_EXTERN PetscErrorCode MatConvert(Mat,MatType,MatReuse,Mat*);
564: PETSC_EXTERN PetscErrorCode MatDuplicate(Mat,MatDuplicateOption,Mat*);


567: PETSC_EXTERN PetscErrorCode MatCopy(Mat,Mat,MatStructure);
568: PETSC_EXTERN PetscErrorCode MatView(Mat,PetscViewer);
569: PETSC_EXTERN PetscErrorCode MatIsSymmetric(Mat,PetscReal,PetscBool *);
570: PETSC_EXTERN PetscErrorCode MatIsStructurallySymmetric(Mat,PetscBool *);
571: PETSC_EXTERN PetscErrorCode MatIsHermitian(Mat,PetscReal,PetscBool *);
572: PETSC_EXTERN PetscErrorCode MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *);
573: PETSC_EXTERN PetscErrorCode MatIsHermitianKnown(Mat,PetscBool *,PetscBool *);
574: PETSC_EXTERN PetscErrorCode MatMissingDiagonal(Mat,PetscBool  *,PetscInt *);
575: PETSC_EXTERN PetscErrorCode MatLoad(Mat, PetscViewer);

577: PETSC_EXTERN PetscErrorCode MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
578: PETSC_EXTERN PetscErrorCode MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool  *);
579: PETSC_EXTERN PetscErrorCode MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
580: PETSC_EXTERN PetscErrorCode MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool  *);

582: /*S
583:      MatInfo - Context of matrix information, used with MatGetInfo()

585:    In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE

587:    Level: intermediate

589: .seealso:  MatGetInfo(), MatInfoType
590: S*/
591: typedef struct {
592:   PetscLogDouble block_size;                         /* block size */
593:   PetscLogDouble nz_allocated,nz_used,nz_unneeded;   /* number of nonzeros */
594:   PetscLogDouble memory;                             /* memory allocated */
595:   PetscLogDouble assemblies;                         /* number of matrix assemblies called */
596:   PetscLogDouble mallocs;                            /* number of mallocs during MatSetValues() */
597:   PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
598:   PetscLogDouble factor_mallocs;                     /* number of mallocs during factorization */
599: } MatInfo;

601: /*E
602:     MatInfoType - Indicates if you want information about the local part of the matrix,
603:      the entire parallel matrix or the maximum over all the local parts.

605:     Level: beginner

607:    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h

609: .seealso: MatGetInfo(), MatInfo
610: E*/
611: typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
612: PETSC_EXTERN PetscErrorCode MatGetInfo(Mat,MatInfoType,MatInfo*);
613: PETSC_EXTERN PetscErrorCode MatGetDiagonal(Mat,Vec);
614: PETSC_EXTERN PetscErrorCode MatGetRowMax(Mat,Vec,PetscInt[]);
615: PETSC_EXTERN PetscErrorCode MatGetRowMin(Mat,Vec,PetscInt[]);
616: PETSC_EXTERN PetscErrorCode MatGetRowMaxAbs(Mat,Vec,PetscInt[]);
617: PETSC_EXTERN PetscErrorCode MatGetRowMinAbs(Mat,Vec,PetscInt[]);
618: PETSC_EXTERN PetscErrorCode MatGetRowSum(Mat,Vec);
619: PETSC_EXTERN PetscErrorCode MatTranspose(Mat,MatReuse,Mat*);
620: PETSC_EXTERN PetscErrorCode MatHermitianTranspose(Mat,MatReuse,Mat*);
621: PETSC_EXTERN PetscErrorCode MatPermute(Mat,IS,IS,Mat*);
622: PETSC_EXTERN PetscErrorCode MatDiagonalScale(Mat,Vec,Vec);
623: PETSC_EXTERN PetscErrorCode MatDiagonalSet(Mat,Vec,InsertMode);

625: PETSC_EXTERN PetscErrorCode MatEqual(Mat,Mat,PetscBool*);
626: PETSC_EXTERN PetscErrorCode MatMultEqual(Mat,Mat,PetscInt,PetscBool*);
627: PETSC_EXTERN PetscErrorCode MatMultAddEqual(Mat,Mat,PetscInt,PetscBool*);
628: PETSC_EXTERN PetscErrorCode MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool*);
629: PETSC_EXTERN PetscErrorCode MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool*);
630: PETSC_EXTERN PetscErrorCode MatMatMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*);
631: PETSC_EXTERN PetscErrorCode MatTransposeMatMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*);
632: PETSC_EXTERN PetscErrorCode MatMatTransposeMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*);
633: PETSC_EXTERN PetscErrorCode MatPtAPMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*);
634: PETSC_EXTERN PetscErrorCode MatRARtMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*);
635: PETSC_EXTERN PetscErrorCode MatIsLinear(Mat,PetscInt,PetscBool*);

637: PETSC_EXTERN PetscErrorCode MatNorm(Mat,NormType,PetscReal*);
638: PETSC_EXTERN PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal*);
639: PETSC_EXTERN PetscErrorCode MatZeroEntries(Mat);
640: PETSC_EXTERN PetscErrorCode MatSetInf(Mat);
641: PETSC_EXTERN PetscErrorCode MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
642: PETSC_EXTERN PetscErrorCode MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec);
643: PETSC_EXTERN PetscErrorCode MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec);
644: PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsStencil(Mat,PetscInt,const MatStencil[],PetscScalar,Vec,Vec);
645: PETSC_EXTERN PetscErrorCode MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
646: PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec);

648: PETSC_EXTERN PetscErrorCode MatGetSize(Mat,PetscInt*,PetscInt*);
649: PETSC_EXTERN PetscErrorCode MatGetLocalSize(Mat,PetscInt*,PetscInt*);
650: PETSC_EXTERN PetscErrorCode MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
651: PETSC_EXTERN PetscErrorCode MatGetOwnershipRanges(Mat,const PetscInt**);
652: PETSC_EXTERN PetscErrorCode MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*);
653: PETSC_EXTERN PetscErrorCode MatGetOwnershipRangesColumn(Mat,const PetscInt**);
654: PETSC_EXTERN PetscErrorCode MatGetOwnershipIS(Mat,IS*,IS*);

656: PETSC_EXTERN PetscErrorCode MatCreateSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
657: PETSC_DEPRECATED_FUNCTION("Use MatCreateSubMatrices() (since version 3.8)") PETSC_STATIC_INLINE PetscErrorCode MatGetSubMatrices(Mat mat,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *submat[]) {return MatCreateSubMatrices(mat,n,irow,icol,scall,submat);}
658: PETSC_EXTERN PetscErrorCode MatCreateSubMatricesMPI(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
659: PETSC_DEPRECATED_FUNCTION("Use MatCreateSubMatricesMPI() (since version 3.8)") PETSC_STATIC_INLINE PetscErrorCode MatGetSubMatricesMPI(Mat mat,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *submat[]) {return MatCreateSubMatricesMPI(mat,n,irow,icol,scall,submat);}
660: PETSC_EXTERN PetscErrorCode MatDestroyMatrices(PetscInt,Mat *[]);
661: PETSC_EXTERN PetscErrorCode MatDestroySubMatrices(PetscInt,Mat *[]);
662: PETSC_EXTERN PetscErrorCode MatCreateSubMatrix(Mat,IS,IS,MatReuse,Mat *);
663: PETSC_DEPRECATED_FUNCTION("Use MatCreateSubMatrix() (since version 3.8)") PETSC_STATIC_INLINE PetscErrorCode MatGetSubMatrix(Mat mat,IS isrow,IS iscol,MatReuse cll,Mat *newmat) {return MatCreateSubMatrix(mat,isrow,iscol,cll,newmat);}
664: PETSC_EXTERN PetscErrorCode MatGetLocalSubMatrix(Mat,IS,IS,Mat*);
665: PETSC_EXTERN PetscErrorCode MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*);
666: PETSC_EXTERN PetscErrorCode MatGetSeqNonzeroStructure(Mat,Mat*);
667: PETSC_EXTERN PetscErrorCode MatDestroySeqNonzeroStructure(Mat*);

669: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJ(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*);
670: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*);
671: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJNumeric(Mat,Mat);
672: PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMat(Mat,MatReuse,Mat*);
673: PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*);
674: PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMatMerge(Mat,MatReuse,IS*,Mat*);
675: PETSC_EXTERN PetscErrorCode MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,Mat*);
676: PETSC_EXTERN PetscErrorCode MatGetGhosts(Mat, PetscInt *,const PetscInt *[]);

678: PETSC_EXTERN PetscErrorCode MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt);
679: PETSC_EXTERN PetscErrorCode MatIncreaseOverlapSplit(Mat mat,PetscInt n,IS is[],PetscInt ov);
680: PETSC_EXTERN PetscErrorCode MatMPIAIJSetUseScalableIncreaseOverlap(Mat,PetscBool);

682: PETSC_EXTERN PetscErrorCode MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);

684: PETSC_EXTERN PetscErrorCode MatMatMatMult(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
685: PETSC_EXTERN PetscErrorCode MatGalerkin(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);

687: PETSC_EXTERN PetscErrorCode MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*);
688: PETSC_EXTERN PetscErrorCode MatRARt(Mat,Mat,MatReuse,PetscReal,Mat*);

690: PETSC_EXTERN PetscErrorCode MatTransposeMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
691: PETSC_EXTERN PetscErrorCode MatMatTransposeMult(Mat,Mat,MatReuse,PetscReal,Mat*);

693: PETSC_EXTERN PetscErrorCode MatAXPY(Mat,PetscScalar,Mat,MatStructure);
694: PETSC_EXTERN PetscErrorCode MatAYPX(Mat,PetscScalar,Mat,MatStructure);

696: PETSC_EXTERN PetscErrorCode MatScale(Mat,PetscScalar);
697: PETSC_EXTERN PetscErrorCode MatShift(Mat,PetscScalar);

699: PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
700: PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*);
701: PETSC_EXTERN PetscErrorCode MatGetLayouts(Mat,PetscLayout*,PetscLayout*);
702: PETSC_EXTERN PetscErrorCode MatSetLayouts(Mat,PetscLayout,PetscLayout);
703: PETSC_EXTERN PetscErrorCode MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
704: PETSC_EXTERN PetscErrorCode MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
705: PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
706: PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
707: PETSC_EXTERN PetscErrorCode MatGetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
708: PETSC_EXTERN PetscErrorCode MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
709: PETSC_EXTERN PetscErrorCode MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);

711: PETSC_EXTERN PetscErrorCode MatStashSetInitialSize(Mat,PetscInt,PetscInt);
712: PETSC_EXTERN PetscErrorCode MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);

714: PETSC_EXTERN PetscErrorCode MatInterpolate(Mat,Vec,Vec);
715: PETSC_EXTERN PetscErrorCode MatInterpolateAdd(Mat,Vec,Vec,Vec);
716: PETSC_EXTERN PetscErrorCode MatRestrict(Mat,Vec,Vec);
717: PETSC_EXTERN PetscErrorCode MatMatInterpolate(Mat,Mat,Mat*);
718: PETSC_EXTERN PetscErrorCode MatMatInterpolateAdd(Mat,Mat,Mat,Mat*);
719: PETSC_EXTERN PetscErrorCode MatMatRestrict(Mat,Mat,Mat*);
720: PETSC_EXTERN PetscErrorCode MatCreateVecs(Mat,Vec*,Vec*);
721: PETSC_DEPRECATED_FUNCTION("Use MatCreateVecs() (since version 3.6)") PETSC_STATIC_INLINE PetscErrorCode MatGetVecs(Mat mat,Vec *x,Vec *y) {return MatCreateVecs(mat,x,y);}
722: PETSC_EXTERN PetscErrorCode MatCreateRedundantMatrix(Mat,PetscInt,MPI_Comm,MatReuse,Mat*);
723: PETSC_EXTERN PetscErrorCode MatGetMultiProcBlock(Mat,MPI_Comm,MatReuse,Mat*);
724: PETSC_EXTERN PetscErrorCode MatFindZeroDiagonals(Mat,IS*);
725: PETSC_EXTERN PetscErrorCode MatFindOffBlockDiagonalEntries(Mat,IS*);
726: PETSC_EXTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);

728: /*MC
729:    MatSetValue - Set a single entry into a matrix.

731:    Not collective

733:    Synopsis:
734: #include <petscmat.h>
735:      PetscErrorCode MatSetValue(Mat m,PetscInt row,PetscInt col,PetscScalar value,InsertMode mode)

737:    Input Parameters:
738: +  m - the matrix
739: .  row - the row location of the entry
740: .  col - the column location of the entry
741: .  value - the value to insert
742: -  mode - either INSERT_VALUES or ADD_VALUES

744:    Notes:
745:    For efficiency one should use MatSetValues() and set several or many
746:    values simultaneously if possible.

748:    Level: beginner

750: .seealso: MatSetValues(), MatSetValueLocal()
751: M*/
752: PETSC_STATIC_INLINE PetscErrorCode MatSetValue(Mat v,PetscInt i,PetscInt j,PetscScalar va,InsertMode mode) {return MatSetValues(v,1,&i,1,&j,&va,mode);}

754: PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);}

756: PETSC_STATIC_INLINE PetscErrorCode MatSetValueLocal(Mat v,PetscInt i,PetscInt j,PetscScalar va,InsertMode mode) {return MatSetValuesLocal(v,1,&i,1,&j,&va,mode);}

758: /*MC
759:    MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per
760:        row in a matrix providing the data that one can use to correctly preallocate the matrix.

762:    Synopsis:
763: #include <petscmat.h>
764:    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)

766:    Collective

768:    Input Parameters:
769: +  comm - the communicator that will share the eventually allocated matrix
770: .  nrows - the number of LOCAL rows in the matrix
771: -  ncols - the number of LOCAL columns in the matrix

773:    Output Parameters:
774: +  dnz - the array that will be passed to the matrix preallocation routines
775: -  onz - the other array passed to the matrix preallocation routines

777:    Level: intermediate

779:    Notes:
780:     See Users-Manual: ch_performance for more details.

782:    Do not malloc or free dnz and onz, that is handled internally by these routines

784:    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().

786: .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
787:           MatPreallocateSymmetricSetLocalBlock()
788: M*/
789: #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
790: do { \
791:   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ncols = (ncols),__rstart,__start,__end; \
792:   _4_PetscCalloc2((size_t)__nrows,&dnz,(size_t)__nrows,&onz);CHKERRQ(_4_ierr); \
793:   __start = 0; __end = __start;                                         \
794:   _4_MPI_Scan(&__ncols,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRMPI(_4_ierr); __start = __end - __ncols;\
795:   _4_MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRMPI(_4_ierr); __rstart = __rstart - __nrows; \
796:   do { } while (0)

798: /*MC
799:    MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
800:        inserted using a local number of the rows and columns

802:    Synopsis:
803: #include <petscmat.h>
804:    PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)

806:    Not Collective

808:    Input Parameters:
809: +  map - the row mapping from local numbering to global numbering
810: .  nrows - the number of rows indicated
811: .  rows - the indices of the rows
812: .  cmap - the column mapping from local to global numbering
813: .  ncols - the number of columns in the matrix
814: .  cols - the columns indicated
815: .  dnz - the array that will be passed to the matrix preallocation routines
816: -  onz - the other array passed to the matrix preallocation routines

818:    Level: intermediate

820:    Notes:
821:     See Users-Manual: ch_performance for more details.

823:    Do not malloc or free dnz and onz, that is handled internally by these routines

825: .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock()
826:           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocalRemoveDups()
827: M*/
828: #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
829: do {\
830:   PetscInt __l;\
831:   _4_ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
832:   _4_ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
833:   for (__l=0;__l<nrows;__l++) {\
834:     _4_MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
835:   }\
836:  } while (0)

838: /*MC
839:    MatPreallocateSetLocalRemoveDups - Indicates the locations (rows and columns) in the matrix where nonzeros will be
840:        inserted using a local number of the rows and columns. This version removes any duplicate columns in cols

842:    Synopsis:
843: #include <petscmat.h>
844:    PetscErrorCode MatPreallocateSetLocalRemoveDups(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)

846:    Not Collective

848:    Input Parameters:
849: +  map - the row mapping from local numbering to global numbering
850: .  nrows - the number of rows indicated
851: .  rows - the indices of the rows (these values are mapped to the global values)
852: .  cmap - the column mapping from local to global numbering
853: .  ncols - the number of columns in the matrix   (this value will be changed if duplicate columns are found)
854: .  cols - the columns indicated (these values are mapped to the global values, they are then sorted and duplicates removed)
855: .  dnz - the array that will be passed to the matrix preallocation routines
856: -  onz - the other array passed to the matrix preallocation routines

858:    Level: intermediate

860:    Notes:
861:     See Users-Manual: ch_performance for more details.

863:    Do not malloc or free dnz and onz, that is handled internally by these routines

865: .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock()
866:           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal()
867: M*/
868: #define MatPreallocateSetLocalRemoveDups(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
869: do {\
870:   PetscInt __l;\
871:   _4_ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
872:   _4_ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
873:   _4_PetscSortRemoveDupsInt(&ncols,cols);CHKERRQ(_4_ierr);\
874:   for (__l=0;__l<nrows;__l++) {\
875:     _4_MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
876:   }\
877:  } while (0)

879: /*MC
880:    MatPreallocateSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be
881:        inserted using a local number of the rows and columns

883:    Synopsis:
884: #include <petscmat.h>
885:    PetscErrorCode MatPreallocateSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)

887:    Not Collective

889:    Input Parameters:
890: +  map - the row mapping from local numbering to global numbering
891: .  nrows - the number of rows indicated
892: .  rows - the indices of the rows
893: .  cmap - the column mapping from local to global numbering
894: .  ncols - the number of columns in the matrix
895: .  cols - the columns indicated
896: .  dnz - the array that will be passed to the matrix preallocation routines
897: -  onz - the other array passed to the matrix preallocation routines

899:    Level: intermediate

901:    Notes:
902:     See Users-Manual: ch_performance for more details.

904:    Do not malloc or free dnz and onz, that is handled internally by these routines

906: .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock()
907:           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock()
908: M*/
909: #define MatPreallocateSetLocalBlock(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
910: do {\
911:   PetscInt __l;\
912:   _4_ISLocalToGlobalMappingApplyBlock(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
913:   _4_ISLocalToGlobalMappingApplyBlock(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
914:   for (__l=0;__l<nrows;__l++) {\
915:     _4_MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
916:   }\
917: } while (0)

919: /*MC
920:    MatPreallocateSymmetricSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be
921:        inserted using a local number of the rows and columns

923:    Synopsis:
924: #include <petscmat.h>
925:    PetscErrorCode MatPreallocateSymmetricSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)

927:    Not Collective

929:    Input Parameters:
930: +  map - the mapping between local numbering and global numbering
931: .  nrows - the number of rows indicated
932: .  rows - the indices of the rows
933: .  ncols - the number of columns in the matrix
934: .  cols - the columns indicated
935: .  dnz - the array that will be passed to the matrix preallocation routines
936: -  onz - the other array passed to the matrix preallocation routines

938:    Level: intermediate

940:    Notes:
941:     See Users-Manual: ch_performance for more details.

943:    Do not malloc or free dnz and onz that is handled internally by these routines

945: .seealso: MatPreallocateFinalize(), MatPreallocateSet()
946:           MatPreallocateInitialize(),  MatPreallocateSetLocal()
947: M*/
948: #define MatPreallocateSymmetricSetLocalBlock(map,nrows,rows,ncols,cols,dnz,onz) 0;\
949: do {\
950:   PetscInt __l;\
951:   _4_ISLocalToGlobalMappingApplyBlock(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
952:   _4_ISLocalToGlobalMappingApplyBlock(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
953:   for (__l=0;__l<nrows;__l++) {\
954:     _4_MatPreallocateSymmetricSetBlock((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
955:   }\
956: } while (0)

958: /*MC
959:    MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
960:        inserted using a local number of the rows and columns

962:    Synopsis:
963: #include <petscmat.h>
964:    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)

966:    Not Collective

968:    Input Parameters:
969: +  row - the row
970: .  ncols - the number of columns in the matrix
971: -  cols - the columns indicated

973:    Output Parameters:
974: +  dnz - the array that will be passed to the matrix preallocation routines
975: -  onz - the other array passed to the matrix preallocation routines

977:    Level: intermediate

979:    Notes:
980:     See Users-Manual: ch_performance for more details.

982:    Do not malloc or free dnz and onz that is handled internally by these routines

984:    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().

986: .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock()
987:           MatPreallocateInitialize(), MatPreallocateSetLocal()
988: M*/
989: #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
990: do { PetscInt __i; \
991:   if (row < __rstart) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Trying to set preallocation for row %D less than first local row %D",row,__rstart);\
992:   if (row >= __rstart+__nrows) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Trying to set preallocation for row %D greater than last local row %D",row,__rstart+__nrows-1);\
993:   for (__i=0; __i<nc; __i++) {\
994:     if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \
995:     else if (dnz[row - __rstart] < __ncols) dnz[row - __rstart]++;\
996:   }\
997: } while (0)

999: /*MC
1000:    MatPreallocateSymmetricSetBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be
1001:        inserted using a local number of the rows and columns

1003:    Synopsis:
1004: #include <petscmat.h>
1005:    PetscErrorCode MatPreallocateSymmetricSetBlock(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)

1007:    Not Collective

1009:    Input Parameters:
1010: +  nrows - the number of rows indicated
1011: .  rows - the indices of the rows
1012: .  ncols - the number of columns in the matrix
1013: .  cols - the columns indicated
1014: .  dnz - the array that will be passed to the matrix preallocation routines
1015: -  onz - the other array passed to the matrix preallocation routines

1017:    Level: intermediate

1019:    Notes:
1020:     See Users-Manual: ch_performance for more details.

1022:    Do not malloc or free dnz and onz that is handled internally by these routines

1024:    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().

1026: .seealso: MatPreallocateFinalize(), MatPreallocateSet(),  MatPreallocateInitialize(),
1027:           MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal()
1028: M*/
1029: #define MatPreallocateSymmetricSetBlock(row,nc,cols,dnz,onz) 0;\
1030: do { PetscInt __i; \
1031:   for (__i=0; __i<nc; __i++) {\
1032:     if (cols[__i] >= __end) onz[row - __rstart]++; \
1033:     else if (cols[__i] >= row && dnz[row - __rstart] < __ncols) dnz[row - __rstart]++;\
1034:   }\
1035: } while (0)

1037: /*MC
1038:    MatPreallocateLocation -  An alternative to MatPreallocateSet() that puts the nonzero locations into the matrix if it exists

1040:    Synopsis:
1041: #include <petscmat.h>
1042:    PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz)

1044:    Not Collective

1046:    Input Parameters:
1047: +  A - matrix
1048: .  row - row where values exist (must be local to this process)
1049: .  ncols - number of columns
1050: .  cols - columns with nonzeros
1051: .  dnz - the array that will be passed to the matrix preallocation routines
1052: -  onz - the other array passed to the matrix preallocation routines

1054:    Level: intermediate

1056:    Notes:
1057:     See Users-Manual: ch_performance for more details.

1059:    Do not malloc or free dnz and onz that is handled internally by these routines

1061:    This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines.

1063: .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
1064:           MatPreallocateSymmetricSetLocalBlock()
1065: M*/
1066: #define MatPreallocateLocation(A,row,ncols,cols,dnz,onz) 0; do {if (A) {MatSetValues(A,1,&row,ncols,cols,NULL,INSERT_VALUES);} else { MatPreallocateSet(row,ncols,cols,dnz,onz);}} while (0)


1069: /*MC
1070:    MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per
1071:        row in a matrix providing the data that one can use to correctly preallocate the matrix.

1073:    Synopsis:
1074: #include <petscmat.h>
1075:    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)

1077:    Collective

1079:    Input Parameters:
1080: +  dnz - the array that was be passed to the matrix preallocation routines
1081: -  onz - the other array passed to the matrix preallocation routines

1083:    Level: intermediate

1085:    Notes:
1086:     See Users-Manual: ch_performance for more details.

1088:    Do not malloc or free dnz and onz that is handled internally by these routines

1090:    This is a MACRO not a function because it closes the { started in MatPreallocateInitialize().

1092: .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
1093:           MatPreallocateSymmetricSetLocalBlock()
1094: M*/
1095: #define MatPreallocateFinalize(dnz,onz) 0;_4_PetscFree2(dnz,onz);CHKERRQ(_4_ierr);} while (0)

1097: /* Routines unique to particular data structures */
1098: PETSC_EXTERN PetscErrorCode MatShellGetContext(Mat,void *);

1100: PETSC_EXTERN PetscErrorCode MatInodeAdjustForInodes(Mat,IS*,IS*);
1101: PETSC_EXTERN PetscErrorCode MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);

1103: PETSC_EXTERN PetscErrorCode MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
1104: PETSC_EXTERN PetscErrorCode MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
1105: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
1106: PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
1107: PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
1108: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJFromTriple(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*,PetscInt,PetscBool);

1110: #define MAT_SKIP_ALLOCATION -4

1112: PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
1113: PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
1114: PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
1115: PETSC_EXTERN PetscErrorCode MatSeqAIJSetTotalPreallocation(Mat,PetscInt);

1117: PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1118: PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1119: PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1120: PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []);
1121: PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
1122: PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
1123: PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
1124: PETSC_EXTERN PetscErrorCode MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
1125: PETSC_EXTERN PetscErrorCode MatMPIAdjToSeq(Mat,Mat*);
1126: PETSC_EXTERN PetscErrorCode MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
1127: PETSC_EXTERN PetscErrorCode MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
1128: PETSC_EXTERN PetscErrorCode MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,const PetscInt*[]);
1129: PETSC_EXTERN PetscErrorCode MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,const PetscInt*[]);
1130: PETSC_EXTERN PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat,Mat*);

1132: PETSC_EXTERN PetscErrorCode MatDenseGetLDA(Mat,PetscInt*);
1133: PETSC_EXTERN PetscErrorCode MatDenseSetLDA(Mat,PetscInt);
1134: PETSC_DEPRECATED_FUNCTION("Use MatDenseSetLDA() (since version 3.14)") PETSC_STATIC_INLINE PetscErrorCode MatSeqDenseSetLDA(Mat A,PetscInt lda) {return MatDenseSetLDA(A,lda);}
1135: PETSC_EXTERN PetscErrorCode MatDenseGetLocalMatrix(Mat,Mat*);

1137: PETSC_EXTERN PetscErrorCode MatBlockMatSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);

1139: PETSC_EXTERN PetscErrorCode MatStoreValues(Mat);
1140: PETSC_EXTERN PetscErrorCode MatRetrieveValues(Mat);

1142: PETSC_EXTERN PetscErrorCode MatFindNonzeroRows(Mat,IS*);
1143: PETSC_EXTERN PetscErrorCode MatFindZeroRows(Mat,IS*);
1144: /*
1145:   These routines are not usually accessed directly, rather solving is
1146:   done through the KSP and PC interfaces.
1147: */

1149: /*J
1150:     MatOrderingType - String with the name of a PETSc matrix ordering

1152:    Level: beginner

1154:    Notes:
1155:       If MATORDERINGEXTERNAL is used then PETSc does not compute an ordering and utilizes one built into the factorization package

1157: .seealso: MatGetOrdering()
1158: J*/
1159: typedef const char* MatOrderingType;
1160: #define MATORDERINGNATURAL        "natural"
1161: #define MATORDERINGND             "nd"
1162: #define MATORDERING1WD            "1wd"
1163: #define MATORDERINGRCM            "rcm"
1164: #define MATORDERINGQMD            "qmd"
1165: #define MATORDERINGROWLENGTH      "rowlength"
1166: #define MATORDERINGWBM            "wbm"
1167: #define MATORDERINGSPECTRAL       "spectral"
1168: #define MATORDERINGAMD            "amd"            /* only works if UMFPACK is installed with PETSc */
1169: #define MATORDERINGNATURAL_OR_ND  "natural_or_nd"  /* special coase used for Cholesky and ICC, allows ND when AIJ matrix is used but Natural when SBAIJ is used */
1170: #define MATORDERINGEXTERNAL       "external"       /* uses an ordering type internal to the factorization package */

1172: PETSC_EXTERN PetscErrorCode MatGetOrdering(Mat,MatOrderingType,IS*,IS*);
1173: PETSC_EXTERN PetscErrorCode MatGetOrderingList(PetscFunctionList*);
1174: PETSC_EXTERN PetscErrorCode MatOrderingRegister(const char[],PetscErrorCode(*)(Mat,MatOrderingType,IS*,IS*));
1175: PETSC_EXTERN PetscFunctionList MatOrderingList;

1177: PETSC_EXTERN PetscErrorCode MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
1178: PETSC_EXTERN PetscErrorCode MatCreateLaplacian(Mat,PetscReal,PetscBool,Mat*);

1180: PETSC_EXTERN PetscErrorCode MatFactorGetPreferredOrdering(Mat,MatFactorType,MatOrderingType*);

1182: /*S
1183:     MatFactorShiftType - Numeric Shift for factorizations

1185:    Level: beginner

1187: .seealso: MatGetFactor()
1188: S*/
1189: typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType;
1190: PETSC_EXTERN const char *const MatFactorShiftTypes[];
1191: PETSC_EXTERN const char *const MatFactorShiftTypesDetail[];

1193: /*S
1194:     MatFactorError - indicates what type of error was generated in a matrix factorization

1196:     Level: beginner

1198:     Developer Notes:
1199:     Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h

1201: .seealso: MatGetFactor()
1202: S*/
1203: typedef enum {MAT_FACTOR_NOERROR,MAT_FACTOR_STRUCT_ZEROPIVOT,MAT_FACTOR_NUMERIC_ZEROPIVOT,MAT_FACTOR_OUTMEMORY,MAT_FACTOR_OTHER} MatFactorError;

1205: PETSC_EXTERN PetscErrorCode MatFactorGetError(Mat,MatFactorError*);
1206: PETSC_EXTERN PetscErrorCode MatFactorClearError(Mat);
1207: PETSC_EXTERN PetscErrorCode MatFactorGetErrorZeroPivot(Mat,PetscReal*,PetscInt*);

1209: /*S
1210:    MatFactorInfo - Data passed into the matrix factorization routines, and information about the resulting factorization

1212:    In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use
1213: $     MatFactorInfo  info(MAT_FACTORINFO_SIZE)

1215:    Notes:
1216:     These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC.

1218:       You can use MatFactorInfoInitialize() to set default values.

1220:    Level: developer

1222: .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
1223:           MatFactorInfoInitialize()

1225: S*/
1226: typedef struct {
1227:   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
1228:   PetscReal     usedt;
1229:   PetscReal     dt;             /* drop tolerance */
1230:   PetscReal     dtcol;          /* tolerance for pivoting */
1231:   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
1232:   PetscReal     fill;           /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */
1233:   PetscReal     levels;         /* ICC/ILU(levels) */
1234:   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
1235:                                    factorization may be faster if do not pivot */
1236:   PetscReal     zeropivot;      /* pivot is called zero if less than this */
1237:   PetscReal     shifttype;      /* type of shift added to matrix factor to prevent zero pivots */
1238:   PetscReal     shiftamount;     /* how large the shift is */
1239: } MatFactorInfo;

1241: PETSC_EXTERN PetscErrorCode MatFactorInfoInitialize(MatFactorInfo*);
1242: PETSC_EXTERN PetscErrorCode MatCholeskyFactor(Mat,IS,const MatFactorInfo*);
1243: PETSC_EXTERN PetscErrorCode MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
1244: PETSC_EXTERN PetscErrorCode MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*);
1245: PETSC_EXTERN PetscErrorCode MatLUFactor(Mat,IS,IS,const MatFactorInfo*);
1246: PETSC_EXTERN PetscErrorCode MatILUFactor(Mat,IS,IS,const MatFactorInfo*);
1247: PETSC_EXTERN PetscErrorCode MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
1248: PETSC_EXTERN PetscErrorCode MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
1249: PETSC_EXTERN PetscErrorCode MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
1250: PETSC_EXTERN PetscErrorCode MatICCFactor(Mat,IS,const MatFactorInfo*);
1251: PETSC_EXTERN PetscErrorCode MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*);
1252: PETSC_EXTERN PetscErrorCode MatQRFactor(Mat,IS,const MatFactorInfo*);
1253: PETSC_EXTERN PetscErrorCode MatQRFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
1254: PETSC_EXTERN PetscErrorCode MatQRFactorNumeric(Mat,Mat,const MatFactorInfo*);
1255: PETSC_EXTERN PetscErrorCode MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
1256: PETSC_EXTERN PetscErrorCode MatSolve(Mat,Vec,Vec);
1257: PETSC_EXTERN PetscErrorCode MatForwardSolve(Mat,Vec,Vec);
1258: PETSC_EXTERN PetscErrorCode MatBackwardSolve(Mat,Vec,Vec);
1259: PETSC_EXTERN PetscErrorCode MatSolveAdd(Mat,Vec,Vec,Vec);
1260: PETSC_EXTERN PetscErrorCode MatSolveTranspose(Mat,Vec,Vec);
1261: PETSC_EXTERN PetscErrorCode MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
1262: PETSC_EXTERN PetscErrorCode MatSolves(Mat,Vecs,Vecs);
1263: PETSC_EXTERN PetscErrorCode MatSetUnfactored(Mat);

1265: typedef enum {MAT_FACTOR_SCHUR_UNFACTORED, MAT_FACTOR_SCHUR_FACTORED, MAT_FACTOR_SCHUR_INVERTED} MatFactorSchurStatus;
1266: PETSC_EXTERN PetscErrorCode MatFactorSetSchurIS(Mat,IS);
1267: PETSC_EXTERN PetscErrorCode MatFactorGetSchurComplement(Mat,Mat*,MatFactorSchurStatus*);
1268: PETSC_EXTERN PetscErrorCode MatFactorRestoreSchurComplement(Mat,Mat*,MatFactorSchurStatus);
1269: PETSC_EXTERN PetscErrorCode MatFactorInvertSchurComplement(Mat);
1270: PETSC_EXTERN PetscErrorCode MatFactorCreateSchurComplement(Mat,Mat*,MatFactorSchurStatus*);
1271: PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplement(Mat,Vec,Vec);
1272: PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplementTranspose(Mat,Vec,Vec);
1273: PETSC_EXTERN PetscErrorCode MatFactorFactorizeSchurComplement(Mat);

1275: /*E
1276:     MatSORType - What type of (S)SOR to perform

1278:     Level: beginner

1280:    May be bitwise ORd together

1282:    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h

1284:    MatSORType may be bitwise ORd together, so do not change the numbers

1286: .seealso: MatSOR()
1287: E*/
1288: typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1289:               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1290:               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
1291:               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
1292: PETSC_EXTERN PetscErrorCode MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);

1294: /*
1295:     These routines are for efficiently computing Jacobians via finite differences.
1296: */

1298: /*S
1299:      MatColoring - Object for managing the coloring of matrices.

1301:    Level: beginner

1303:     Notes:
1304:        Coloring of matrices can be computed directly from the sparse matrix nonzero structure via the MatColoring object or from the mesh from which the
1305:        matrix comes from via DMCreateColoring(). In general using the mesh produces a more optimal coloring (fewer colors).

1307:        Once a coloring is available MatFDColoringCreate() creates an object that can be used to efficiently compute Jacobians using that coloring. This
1308:        same object can also be used to efficiently convert data created by Automatic Differentation tools to PETSc sparse matrices.

1310: .seealso:  MatFDColoringCreate(), MatColoringWeightType, ISColoring, MatFDColoring, DMCreateColoring(), MatColoringCreate(), MatOrdering, MatPartitioning
1311: S*/
1312: typedef struct _p_MatColoring* MatColoring;

1314: /*J
1315:     MatColoringType - String with the name of a PETSc matrix coloring

1317:    Level: beginner

1319: .seealso: MatColoringSetType(), MatColoring
1320: J*/
1321: typedef const  char*           MatColoringType;
1322: #define MATCOLORINGJP      "jp"
1323: #define MATCOLORINGPOWER   "power"
1324: #define MATCOLORINGNATURAL "natural"
1325: #define MATCOLORINGSL      "sl"
1326: #define MATCOLORINGLF      "lf"
1327: #define MATCOLORINGID      "id"
1328: #define MATCOLORINGGREEDY  "greedy"

1330: /*E
1331:    MatColoringWeightType - Type of weight scheme

1333:     Not Collective

1335: +   MAT_COLORING_RANDOM  - Random weights
1336: .   MAT_COLORING_LEXICAL - Lexical weighting based upon global numbering.
1337: -   MAT_COLORING_LF      - Last-first weighting.

1339:     Level: intermediate

1341:    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h

1343: .seealso: MatColoring, MatColoringCreate()
1344: E*/
1345: typedef enum {MAT_COLORING_WEIGHT_RANDOM,MAT_COLORING_WEIGHT_LEXICAL,MAT_COLORING_WEIGHT_LF,MAT_COLORING_WEIGHT_SL} MatColoringWeightType;

1347: PETSC_EXTERN PetscErrorCode MatColoringCreate(Mat,MatColoring*);
1348: PETSC_EXTERN PetscErrorCode MatColoringGetDegrees(Mat,PetscInt,PetscInt*);
1349: PETSC_EXTERN PetscErrorCode MatColoringDestroy(MatColoring*);
1350: PETSC_EXTERN PetscErrorCode MatColoringView(MatColoring,PetscViewer);
1351: PETSC_EXTERN PetscErrorCode MatColoringSetType(MatColoring,MatColoringType);
1352: PETSC_EXTERN PetscErrorCode MatColoringSetFromOptions(MatColoring);
1353: PETSC_EXTERN PetscErrorCode MatColoringSetDistance(MatColoring,PetscInt);
1354: PETSC_EXTERN PetscErrorCode MatColoringGetDistance(MatColoring,PetscInt*);
1355: PETSC_EXTERN PetscErrorCode MatColoringSetMaxColors(MatColoring,PetscInt);
1356: PETSC_EXTERN PetscErrorCode MatColoringGetMaxColors(MatColoring,PetscInt*);
1357: PETSC_EXTERN PetscErrorCode MatColoringApply(MatColoring,ISColoring*);
1358: PETSC_EXTERN PetscErrorCode MatColoringRegister(const char[],PetscErrorCode(*)(MatColoring));
1359: PETSC_EXTERN PetscErrorCode MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
1360: PETSC_EXTERN PetscErrorCode MatColoringSetWeightType(MatColoring,MatColoringWeightType);
1361: PETSC_EXTERN PetscErrorCode MatColoringSetWeights(MatColoring,PetscReal*,PetscInt*);
1362: PETSC_EXTERN PetscErrorCode MatColoringCreateWeights(MatColoring,PetscReal **,PetscInt **lperm);
1363: PETSC_EXTERN PetscErrorCode MatColoringTest(MatColoring,ISColoring);
1364: PETSC_DEPRECATED_FUNCTION("Use MatColoringTest() (since version 3.10)") PETSC_STATIC_INLINE PetscErrorCode MatColoringTestValid(MatColoring matcoloring,ISColoring iscoloring) {return MatColoringTest(matcoloring,iscoloring);}
1365: PETSC_EXTERN PetscErrorCode MatISColoringTest(Mat,ISColoring);

1367: /*S
1368:      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1369:         and coloring

1371:    Level: beginner

1373:    Notes:
1374:       This object is creating utilizing a coloring provided by the MatColoring object or DMCreateColoring()

1376: .seealso:  MatFDColoringCreate(), MatColoring, DMCreateColoring()
1377: S*/
1378: typedef struct _p_MatFDColoring* MatFDColoring;

1380: PETSC_EXTERN PetscErrorCode MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
1381: PETSC_EXTERN PetscErrorCode MatFDColoringDestroy(MatFDColoring*);
1382: PETSC_EXTERN PetscErrorCode MatFDColoringView(MatFDColoring,PetscViewer);
1383: PETSC_EXTERN PetscErrorCode MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
1384: PETSC_EXTERN PetscErrorCode MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**);
1385: PETSC_EXTERN PetscErrorCode MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
1386: PETSC_EXTERN PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring);
1387: PETSC_EXTERN PetscErrorCode MatFDColoringApply(Mat,MatFDColoring,Vec,void *);
1388: PETSC_EXTERN PetscErrorCode MatFDColoringSetF(MatFDColoring,Vec);
1389: PETSC_EXTERN PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,const PetscInt*[]);
1390: PETSC_EXTERN PetscErrorCode MatFDColoringSetUp(Mat,ISColoring,MatFDColoring);
1391: PETSC_EXTERN PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring,PetscInt,PetscInt);
1392: PETSC_EXTERN PetscErrorCode MatFDColoringSetValues(Mat,MatFDColoring,const PetscScalar*);

1394: /*S
1395:      MatTransposeColoring - Object for computing a sparse matrix product C=A*B^T via coloring

1397:    Level: beginner

1399: .seealso:  MatTransposeColoringCreate()
1400: S*/
1401: typedef struct _p_MatTransposeColoring* MatTransposeColoring;

1403: PETSC_EXTERN PetscErrorCode MatTransposeColoringCreate(Mat,ISColoring,MatTransposeColoring *);
1404: PETSC_EXTERN PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring,Mat,Mat);
1405: PETSC_EXTERN PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring,Mat,Mat);
1406: PETSC_EXTERN PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring*);

1408: /*
1409:     These routines are for partitioning matrices: currently used only
1410:   for adjacency matrix, MatCreateMPIAdj().
1411: */

1413: /*S
1414:      MatPartitioning - Object for managing the partitioning of a matrix or graph

1416:    Level: beginner

1418:    Notes:
1419:      There is also a PetscPartitioner object that provides the same functionality. It can utilize the MatPartitioning operations
1420:      via PetscPartitionerSetType(p,PETSCPARTITIONERMATPARTITIONING)

1422:    Developers Note:
1423:      It is an extra maintainance and documentation cost to have two objects with the same functionality.

1425: .seealso:  MatPartitioningCreate(), MatPartitioningType, MatColoring, MatGetOrdering()
1426: S*/
1427: typedef struct _p_MatPartitioning* MatPartitioning;

1429: /*J
1430:     MatPartitioningType - String with the name of a PETSc matrix partitioning

1432:    Level: beginner
1433: dm
1434: .seealso: MatPartitioningCreate(), MatPartitioning
1435: J*/
1436: typedef const char* MatPartitioningType;
1437: #define MATPARTITIONINGCURRENT  "current"
1438: #define MATPARTITIONINGAVERAGE   "average"
1439: #define MATPARTITIONINGSQUARE   "square"
1440: #define MATPARTITIONINGPARMETIS "parmetis"
1441: #define MATPARTITIONINGCHACO    "chaco"
1442: #define MATPARTITIONINGPARTY    "party"
1443: #define MATPARTITIONINGPTSCOTCH "ptscotch"
1444: #define MATPARTITIONINGHIERARCH  "hierarch"

1446: PETSC_EXTERN PetscErrorCode MatPartitioningCreate(MPI_Comm,MatPartitioning*);
1447: PETSC_EXTERN PetscErrorCode MatPartitioningSetType(MatPartitioning,MatPartitioningType);
1448: PETSC_EXTERN PetscErrorCode MatPartitioningSetNParts(MatPartitioning,PetscInt);
1449: PETSC_EXTERN PetscErrorCode MatPartitioningSetAdjacency(MatPartitioning,Mat);
1450: PETSC_EXTERN PetscErrorCode MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
1451: PETSC_EXTERN PetscErrorCode MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
1452: PETSC_EXTERN PetscErrorCode MatPartitioningSetUseEdgeWeights(MatPartitioning,PetscBool);
1453: PETSC_EXTERN PetscErrorCode MatPartitioningGetUseEdgeWeights(MatPartitioning,PetscBool*);
1454: PETSC_EXTERN PetscErrorCode MatPartitioningApply(MatPartitioning,IS*);
1455: PETSC_EXTERN PetscErrorCode MatPartitioningImprove(MatPartitioning,IS*);
1456: PETSC_EXTERN PetscErrorCode MatPartitioningViewImbalance(MatPartitioning,IS);
1457: PETSC_EXTERN PetscErrorCode MatPartitioningApplyND(MatPartitioning,IS*);
1458: PETSC_EXTERN PetscErrorCode MatPartitioningDestroy(MatPartitioning*);
1459: PETSC_EXTERN PetscErrorCode MatPartitioningRegister(const char[],PetscErrorCode (*)(MatPartitioning));
1460: PETSC_EXTERN PetscErrorCode MatPartitioningView(MatPartitioning,PetscViewer);
1461: PETSC_EXTERN PetscErrorCode MatPartitioningViewFromOptions(MatPartitioning,PetscObject,const char[]);
1462: PETSC_EXTERN PetscErrorCode MatPartitioningSetFromOptions(MatPartitioning);
1463: PETSC_EXTERN PetscErrorCode MatPartitioningGetType(MatPartitioning,MatPartitioningType*);

1465: PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetRepartition(MatPartitioning part);
1466: PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
1467: PETSC_EXTERN PetscErrorCode MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *);

1469: typedef enum { MP_CHACO_MULTILEVEL=1,MP_CHACO_SPECTRAL=2,MP_CHACO_LINEAR=4,MP_CHACO_RANDOM=5,MP_CHACO_SCATTERED=6 } MPChacoGlobalType;
1470: PETSC_EXTERN const char *const MPChacoGlobalTypes[];
1471: typedef enum { MP_CHACO_KERNIGHAN=1,MP_CHACO_NONE=2 } MPChacoLocalType;
1472: PETSC_EXTERN const char *const MPChacoLocalTypes[];
1473: typedef enum { MP_CHACO_LANCZOS=0,MP_CHACO_RQI=1 } MPChacoEigenType;
1474: PETSC_EXTERN const char *const MPChacoEigenTypes[];

1476: PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetGlobal(MatPartitioning,MPChacoGlobalType);
1477: PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetGlobal(MatPartitioning,MPChacoGlobalType*);
1478: PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetLocal(MatPartitioning,MPChacoLocalType);
1479: PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetLocal(MatPartitioning,MPChacoLocalType*);
1480: PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
1481: PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
1482: PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenSolver(MatPartitioning,MPChacoEigenType*);
1483: PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenTol(MatPartitioning,PetscReal);
1484: PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenTol(MatPartitioning,PetscReal*);
1485: PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenNumber(MatPartitioning,PetscInt);
1486: PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenNumber(MatPartitioning,PetscInt*);

1488: #define MP_PARTY_OPT "opt"
1489: #define MP_PARTY_LIN "lin"
1490: #define MP_PARTY_SCA "sca"
1491: #define MP_PARTY_RAN "ran"
1492: #define MP_PARTY_GBF "gbf"
1493: #define MP_PARTY_GCF "gcf"
1494: #define MP_PARTY_BUB "bub"
1495: #define MP_PARTY_DEF "def"
1496: PETSC_EXTERN PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning,const char*);
1497: #define MP_PARTY_HELPFUL_SETS "hs"
1498: #define MP_PARTY_KERNIGHAN_LIN "kl"
1499: #define MP_PARTY_NONE "no"
1500: PETSC_EXTERN PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning,const char*);
1501: PETSC_EXTERN PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
1502: PETSC_EXTERN PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning,PetscBool);
1503: PETSC_EXTERN PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool);

1505: typedef enum { MP_PTSCOTCH_DEFAULT,MP_PTSCOTCH_QUALITY,MP_PTSCOTCH_SPEED,MP_PTSCOTCH_BALANCE,MP_PTSCOTCH_SAFETY,MP_PTSCOTCH_SCALABILITY } MPPTScotchStrategyType;
1506: PETSC_EXTERN const char *const MPPTScotchStrategyTypes[];

1508: PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning,PetscReal);
1509: PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning,PetscReal*);
1510: PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning,MPPTScotchStrategyType);
1511: PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning,MPPTScotchStrategyType*);

1513: /*
1514:  * hierarchical partitioning
1515:  */
1516: PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetFineparts(MatPartitioning,IS*);
1517: PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetCoarseparts(MatPartitioning,IS*);
1518: PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNcoarseparts(MatPartitioning,PetscInt);
1519: PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNfineparts(MatPartitioning, PetscInt);

1521: PETSC_EXTERN PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*);
1522: PETSC_EXTERN PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*);

1524: /*
1525:     If you add entries here you must also add them to include/petscmat.h
1526:     and src/mat/f90-mod/petscmat.h
1527: */
1528: typedef enum { MATOP_SET_VALUES=0,
1529:                MATOP_GET_ROW=1,
1530:                MATOP_RESTORE_ROW=2,
1531:                MATOP_MULT=3,
1532:                MATOP_MULT_ADD=4,
1533:                MATOP_MULT_TRANSPOSE=5,
1534:                MATOP_MULT_TRANSPOSE_ADD=6,
1535:                MATOP_SOLVE=7,
1536:                MATOP_SOLVE_ADD=8,
1537:                MATOP_SOLVE_TRANSPOSE=9,
1538:                MATOP_SOLVE_TRANSPOSE_ADD=10,
1539:                MATOP_LUFACTOR=11,
1540:                MATOP_CHOLESKYFACTOR=12,
1541:                MATOP_SOR=13,
1542:                MATOP_TRANSPOSE=14,
1543:                MATOP_GETINFO=15,
1544:                MATOP_EQUAL=16,
1545:                MATOP_GET_DIAGONAL=17,
1546:                MATOP_DIAGONAL_SCALE=18,
1547:                MATOP_NORM=19,
1548:                MATOP_ASSEMBLY_BEGIN=20,
1549:                MATOP_ASSEMBLY_END=21,
1550:                MATOP_SET_OPTION=22,
1551:                MATOP_ZERO_ENTRIES=23,
1552:                MATOP_ZERO_ROWS=24,
1553:                MATOP_LUFACTOR_SYMBOLIC=25,
1554:                MATOP_LUFACTOR_NUMERIC=26,
1555:                MATOP_CHOLESKY_FACTOR_SYMBOLIC=27,
1556:                MATOP_CHOLESKY_FACTOR_NUMERIC=28,
1557:                MATOP_SETUP_PREALLOCATION=29,
1558:                MATOP_ILUFACTOR_SYMBOLIC=30,
1559:                MATOP_ICCFACTOR_SYMBOLIC=31,
1560:                MATOP_GET_DIAGONAL_BLOCK=32,
1561:                MATOP_FREE_INTER_STRUCT=33,
1562:                MATOP_DUPLICATE=34,
1563:                MATOP_FORWARD_SOLVE=35,
1564:                MATOP_BACKWARD_SOLVE=36,
1565:                MATOP_ILUFACTOR=37,
1566:                MATOP_ICCFACTOR=38,
1567:                MATOP_AXPY=39,
1568:                MATOP_CREATE_SUBMATRICES=40,
1569:                MATOP_INCREASE_OVERLAP=41,
1570:                MATOP_GET_VALUES=42,
1571:                MATOP_COPY=43,
1572:                MATOP_GET_ROW_MAX=44,
1573:                MATOP_SCALE=45,
1574:                MATOP_SHIFT=46,
1575:                MATOP_DIAGONAL_SET=47,
1576:                MATOP_ZERO_ROWS_COLUMNS=48,
1577:                MATOP_SET_RANDOM=49,
1578:                MATOP_GET_ROW_IJ=50,
1579:                MATOP_RESTORE_ROW_IJ=51,
1580:                MATOP_GET_COLUMN_IJ=52,
1581:                MATOP_RESTORE_COLUMN_IJ=53,
1582:                MATOP_FDCOLORING_CREATE=54,
1583:                MATOP_COLORING_PATCH=55,
1584:                MATOP_SET_UNFACTORED=56,
1585:                MATOP_PERMUTE=57,
1586:                MATOP_SET_VALUES_BLOCKED=58,
1587:                MATOP_CREATE_SUBMATRIX=59,
1588:                MATOP_DESTROY=60,
1589:                MATOP_VIEW=61,
1590:                MATOP_CONVERT_FROM=62,
1591:                MATOP_MATMAT_MULT=63,
1592:                MATOP_MATMAT_MULT_SYMBOLIC=64,
1593:                MATOP_MATMAT_MULT_NUMERIC=65,
1594:                MATOP_SET_LOCAL_TO_GLOBAL_MAP=66,
1595:                MATOP_SET_VALUES_LOCAL=67,
1596:                MATOP_ZERO_ROWS_LOCAL=68,
1597:                MATOP_GET_ROW_MAX_ABS=69,
1598:                MATOP_GET_ROW_MIN_ABS=70,
1599:                MATOP_CONVERT=71,
1600:                MATOP_SET_COLORING=72,
1601:                /* MATOP_PLACEHOLDER_73=73, */
1602:                MATOP_SET_VALUES_ADIFOR=74,
1603:                MATOP_FD_COLORING_APPLY=75,
1604:                MATOP_SET_FROM_OPTIONS=76,
1605:                MATOP_MULT_CONSTRAINED=77,
1606:                MATOP_MULT_TRANSPOSE_CONSTRAIN=78,
1607:                MATOP_FIND_ZERO_DIAGONALS=79,
1608:                MATOP_MULT_MULTIPLE=80,
1609:                MATOP_SOLVE_MULTIPLE=81,
1610:                MATOP_GET_INERTIA=82,
1611:                MATOP_LOAD=83,
1612:                MATOP_IS_SYMMETRIC=84,
1613:                MATOP_IS_HERMITIAN=85,
1614:                MATOP_IS_STRUCTURALLY_SYMMETRIC=86,
1615:                MATOP_SET_VALUES_BLOCKEDLOCAL=87,
1616:                MATOP_CREATE_VECS=88,
1617:                MATOP_MAT_MULT=89,
1618:                MATOP_MAT_MULT_SYMBOLIC=90,
1619:                MATOP_MAT_MULT_NUMERIC=91,
1620:                MATOP_PTAP=92,
1621:                MATOP_PTAP_SYMBOLIC=93,
1622:                MATOP_PTAP_NUMERIC=94,
1623:                MATOP_MAT_TRANSPOSE_MULT=95,
1624:                MATOP_MAT_TRANSPOSE_MULT_SYMBO=96,
1625:                MATOP_MAT_TRANSPOSE_MULT_NUMER=97,
1626:                /* MATOP_PLACEHOLDER_98=98, */
1627:                MATOP_PRODUCTSETFROMOPTIONS=99,
1628:                MATOP_PRODUCTSYMBOLIC=100,
1629:                MATOP_PRODUCTNUMERIC=101,
1630:                MATOP_CONJUGATE=102,
1631:                /* MATOP_PLACEHOLDER_103=103, */
1632:                MATOP_SET_VALUES_ROW=104,
1633:                MATOP_REAL_PART=105,
1634:                MATOP_IMAGINARY_PART=106,
1635:                MATOP_GET_ROW_UPPER_TRIANGULAR=107,
1636:                MATOP_RESTORE_ROW_UPPER_TRIANG=108,
1637:                MATOP_MAT_SOLVE=109,
1638:                MATOP_MAT_SOLVE_TRANSPOSE=110,
1639:                MATOP_GET_ROW_MIN=111,
1640:                MATOP_GET_COLUMN_VECTOR=112,
1641:                MATOP_MISSING_DIAGONAL=113,
1642:                MATOP_GET_SEQ_NONZERO_STRUCTUR=114,
1643:                MATOP_CREATE=115,
1644:                MATOP_GET_GHOSTS=116,
1645:                MATOP_GET_LOCAL_SUB_MATRIX=117,
1646:                MATOP_RESTORE_LOCALSUB_MATRIX=118,
1647:                MATOP_MULT_DIAGONAL_BLOCK=119,
1648:                MATOP_HERMITIAN_TRANSPOSE=120,
1649:                MATOP_MULT_HERMITIAN_TRANSPOSE=121,
1650:                MATOP_MULT_HERMITIAN_TRANS_ADD=122,
1651:                MATOP_GET_MULTI_PROC_BLOCK=123,
1652:                MATOP_FIND_NONZERO_ROWS=124,
1653:                MATOP_GET_COLUMN_NORMS=125,
1654:                MATOP_INVERT_BLOCK_DIAGONAL=126,
1655:                /* MATOP_PLACEHOLDER_127=127, */
1656:                MATOP_CREATE_SUB_MATRICES_MPI=128,
1657:                MATOP_SET_VALUES_BATCH=129,
1658:                MATOP_TRANSPOSE_MAT_MULT=130,
1659:                MATOP_TRANSPOSE_MAT_MULT_SYMBO=131,
1660:                MATOP_TRANSPOSE_MAT_MULT_NUMER=132,
1661:                MATOP_TRANSPOSE_COLORING_CREAT=133,
1662:                MATOP_TRANS_COLORING_APPLY_SPT=134,
1663:                MATOP_TRANS_COLORING_APPLY_DEN=135,
1664:                MATOP_RART=136,
1665:                MATOP_RART_SYMBOLIC=137,
1666:                MATOP_RART_NUMERIC=138,
1667:                MATOP_SET_BLOCK_SIZES=139,
1668:                MATOP_AYPX=140,
1669:                MATOP_RESIDUAL=141,
1670:                MATOP_FDCOLORING_SETUP=142,
1671:                MATOP_MPICONCATENATESEQ=144,
1672:                MATOP_DESTROYSUBMATRICES=145,
1673:                MATOP_TRANSPOSE_SOLVE=146,
1674:                MATOP_GET_VALUES_LOCAL=147
1675:              } MatOperation;
1676: PETSC_EXTERN PetscErrorCode MatSetOperation(Mat,MatOperation,void(*)(void));
1677: PETSC_EXTERN PetscErrorCode MatGetOperation(Mat,MatOperation,void(**)(void));
1678: PETSC_EXTERN PetscErrorCode MatHasOperation(Mat,MatOperation,PetscBool*);
1679: PETSC_EXTERN PetscErrorCode MatHasCongruentLayouts(Mat,PetscBool*);
1680: PETSC_DEPRECATED_FUNCTION("Use MatProductClear() (since version 3.14)") PETSC_STATIC_INLINE PetscErrorCode MatFreeIntermediateDataStructures(Mat A) {return MatProductClear(A);}
1681: PETSC_EXTERN PetscErrorCode MatShellSetOperation(Mat,MatOperation,void(*)(void));
1682: PETSC_EXTERN PetscErrorCode MatShellGetOperation(Mat,MatOperation,void(**)(void));
1683: PETSC_EXTERN PetscErrorCode MatShellSetContext(Mat,void*);
1684: PETSC_EXTERN PetscErrorCode MatShellSetVecType(Mat,VecType);
1685: PETSC_EXTERN PetscErrorCode MatShellTestMult(Mat,PetscErrorCode (*)(void*,Vec,Vec),Vec,void*,PetscBool*);
1686: PETSC_EXTERN PetscErrorCode MatShellTestMultTranspose(Mat,PetscErrorCode (*)(void*,Vec,Vec),Vec,void*,PetscBool*);
1687: PETSC_EXTERN PetscErrorCode MatShellSetManageScalingShifts(Mat);
1688: PETSC_EXTERN PetscErrorCode MatShellSetMatProductOperation(Mat,MatProductType,PetscErrorCode (*)(Mat,Mat,Mat,void**),PetscErrorCode (*)(Mat,Mat,Mat,void*),PetscErrorCode (*)(void*),MatType,MatType);
1689: PETSC_EXTERN PetscErrorCode MatIsShell(Mat,PetscBool*);

1691: /*
1692:    Codes for matrices stored on disk. By default they are
1693:    stored in a universal format. By changing the format with
1694:    PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will
1695:    be stored in a way natural for the matrix, for example dense matrices
1696:    would be stored as dense. Matrices stored this way may only be
1697:    read into matrices of the same type.
1698: */
1699: #define MATRIX_BINARY_FORMAT_DENSE -1

1701: PETSC_EXTERN PetscErrorCode MatMPIBAIJSetHashTableFactor(Mat,PetscReal);

1703: PETSC_EXTERN PetscErrorCode MatISSetLocalMatType(Mat,MatType);
1704: PETSC_EXTERN PetscErrorCode MatISSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1705: PETSC_EXTERN PetscErrorCode MatISStoreL2L(Mat,PetscBool);
1706: PETSC_EXTERN PetscErrorCode MatISFixLocalEmpty(Mat,PetscBool);
1707: PETSC_EXTERN PetscErrorCode MatISGetLocalMat(Mat,Mat*);
1708: PETSC_EXTERN PetscErrorCode MatISRestoreLocalMat(Mat,Mat*);
1709: PETSC_EXTERN PetscErrorCode MatISSetLocalMat(Mat,Mat);
1710: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use the MatConvert() interface (since version 3.10)") PetscErrorCode MatISGetMPIXAIJ(Mat,MatReuse,Mat*);

1712: /*S
1713:      MatNullSpace - Object that removes a null space from a vector, i.e.
1714:          orthogonalizes the vector to a subsapce

1716:    Level: advanced

1718:   Users manual sections:
1719: .   sec_singular

1721: .seealso:  MatNullSpaceCreate()
1722: S*/
1723: typedef struct _p_MatNullSpace* MatNullSpace;

1725: PETSC_EXTERN PetscErrorCode MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*);
1726: PETSC_EXTERN PetscErrorCode MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*);
1727: PETSC_EXTERN PetscErrorCode MatNullSpaceDestroy(MatNullSpace*);
1728: PETSC_EXTERN PetscErrorCode MatNullSpaceRemove(MatNullSpace,Vec);
1729: PETSC_EXTERN PetscErrorCode MatGetNullSpace(Mat, MatNullSpace *);
1730: PETSC_EXTERN PetscErrorCode MatGetTransposeNullSpace(Mat, MatNullSpace *);
1731: PETSC_EXTERN PetscErrorCode MatSetTransposeNullSpace(Mat,MatNullSpace);
1732: PETSC_EXTERN PetscErrorCode MatSetNullSpace(Mat,MatNullSpace);
1733: PETSC_EXTERN PetscErrorCode MatSetNearNullSpace(Mat,MatNullSpace);
1734: PETSC_EXTERN PetscErrorCode MatGetNearNullSpace(Mat,MatNullSpace*);
1735: PETSC_EXTERN PetscErrorCode MatNullSpaceTest(MatNullSpace,Mat,PetscBool  *);
1736: PETSC_EXTERN PetscErrorCode MatNullSpaceView(MatNullSpace,PetscViewer);
1737: PETSC_EXTERN PetscErrorCode MatNullSpaceGetVecs(MatNullSpace,PetscBool*,PetscInt*,const Vec**);
1738: PETSC_EXTERN PetscErrorCode MatNullSpaceCreateRigidBody(Vec,MatNullSpace*);

1740: PETSC_EXTERN PetscErrorCode MatReorderingSeqSBAIJ(Mat,IS);
1741: PETSC_EXTERN PetscErrorCode MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
1742: PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);

1744: PETSC_EXTERN PetscErrorCode MatCreateMAIJ(Mat,PetscInt,Mat*);
1745: PETSC_EXTERN PetscErrorCode MatMAIJRedimension(Mat,PetscInt,Mat*);
1746: PETSC_EXTERN PetscErrorCode MatMAIJGetAIJ(Mat,Mat*);

1748: PETSC_EXTERN PetscErrorCode MatComputeOperator(Mat,MatType,Mat*);
1749: PETSC_EXTERN PetscErrorCode MatComputeOperatorTranspose(Mat,MatType,Mat*);

1751: PETSC_DEPRECATED_FUNCTION("Use MatComputeOperator() (since version 3.12)") PETSC_STATIC_INLINE PetscErrorCode MatComputeExplicitOperator(Mat A,Mat* B) { return MatComputeOperator(A,NULL,B); }
1752: PETSC_DEPRECATED_FUNCTION("Use MatComputeOperatorTranspose() (since version 3.12)") PETSC_STATIC_INLINE PetscErrorCode MatComputeExplicitOperatorTranspose(Mat A,Mat* B) { return MatComputeOperatorTranspose(A,NULL,B); }

1754: PETSC_EXTERN PetscErrorCode MatCreateKAIJ(Mat,PetscInt,PetscInt,const PetscScalar[],const PetscScalar[],Mat*);
1755: PETSC_EXTERN PetscErrorCode MatKAIJGetAIJ(Mat,Mat*);
1756: PETSC_EXTERN PetscErrorCode MatKAIJGetS(Mat,PetscInt*,PetscInt*,PetscScalar**);
1757: PETSC_EXTERN PetscErrorCode MatKAIJGetSRead(Mat,PetscInt*,PetscInt*,const PetscScalar**);
1758: PETSC_EXTERN PetscErrorCode MatKAIJRestoreS(Mat,PetscScalar**);
1759: PETSC_EXTERN PetscErrorCode MatKAIJRestoreSRead(Mat,const PetscScalar**);
1760: PETSC_EXTERN PetscErrorCode MatKAIJGetT(Mat,PetscInt*,PetscInt*,PetscScalar**);
1761: PETSC_EXTERN PetscErrorCode MatKAIJGetTRead(Mat,PetscInt*,PetscInt*,const PetscScalar**);
1762: PETSC_EXTERN PetscErrorCode MatKAIJRestoreT(Mat,PetscScalar**);
1763: PETSC_EXTERN PetscErrorCode MatKAIJRestoreTRead(Mat,const PetscScalar**);
1764: PETSC_EXTERN PetscErrorCode MatKAIJSetAIJ(Mat,Mat);
1765: PETSC_EXTERN PetscErrorCode MatKAIJSetS(Mat,PetscInt,PetscInt,const PetscScalar[]);
1766: PETSC_EXTERN PetscErrorCode MatKAIJSetT(Mat,PetscInt,PetscInt,const PetscScalar[]);
1767: PETSC_EXTERN PetscErrorCode MatKAIJGetScaledIdentity(Mat,PetscBool*);

1769: PETSC_EXTERN PetscErrorCode MatDiagonalScaleLocal(Mat,Vec);

1771: PETSC_EXTERN PetscErrorCode MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
1772: PETSC_EXTERN PetscErrorCode MatMFFDSetBase(Mat,Vec,Vec);
1773: PETSC_EXTERN PetscErrorCode MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*);
1774: PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*));
1775: PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec));
1776: PETSC_EXTERN PetscErrorCode MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt);
1777: PETSC_EXTERN PetscErrorCode MatMFFDResetHHistory(Mat);
1778: PETSC_EXTERN PetscErrorCode MatMFFDSetFunctionError(Mat,PetscReal);
1779: PETSC_EXTERN PetscErrorCode MatMFFDSetPeriod(Mat,PetscInt);
1780: PETSC_EXTERN PetscErrorCode MatMFFDGetH(Mat,PetscScalar *);
1781: PETSC_EXTERN PetscErrorCode MatMFFDSetOptionsPrefix(Mat,const char[]);
1782: PETSC_EXTERN PetscErrorCode MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*);
1783: PETSC_EXTERN PetscErrorCode MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*);

1785: /*S
1786:     MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free
1787:               Jacobian vector products

1789:     Notes:
1790:     MATMFFD is a specific MatType which uses the MatMFFD data structure

1792:            MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure

1794:     Level: developer

1796: .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister()
1797: S*/
1798: typedef struct _p_MatMFFD* MatMFFD;

1800: /*J
1801:     MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function

1803:    Level: beginner

1805: .seealso: MatMFFDSetType(), MatMFFDRegister()
1806: J*/
1807: typedef const char* MatMFFDType;
1808: #define MATMFFD_DS  "ds"
1809: #define MATMFFD_WP  "wp"

1811: PETSC_EXTERN PetscErrorCode MatMFFDSetType(Mat,MatMFFDType);
1812: PETSC_EXTERN PetscErrorCode MatMFFDRegister(const char[],PetscErrorCode (*)(MatMFFD));

1814: PETSC_EXTERN PetscErrorCode MatMFFDDSSetUmin(Mat,PetscReal);
1815: PETSC_EXTERN PetscErrorCode MatMFFDWPSetComputeNormU(Mat,PetscBool);

1817: PETSC_EXTERN PetscErrorCode MatFDColoringSetType(MatFDColoring,MatMFFDType);

1819: PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
1820: PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);

1822: #ifdef PETSC_HAVE_HARA
1823: PETSC_EXTERN_TYPEDEF typedef PetscScalar (*MatHaraKernel)(PetscInt,PetscReal[],PetscReal[],void*);
1824: PETSC_EXTERN PetscErrorCode MatCreateHaraFromKernel(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscReal[],MatHaraKernel,void*,PetscReal,PetscInt,PetscInt,Mat*);
1825: PETSC_EXTERN PetscErrorCode MatCreateHaraFromMat(Mat,PetscInt,const PetscReal[],PetscReal,PetscInt,PetscInt,PetscInt,PetscReal,Mat*);
1826: PETSC_EXTERN PetscErrorCode MatHaraSetSamplingMat(Mat,Mat,PetscInt,PetscReal);
1827: #endif

1829: /*
1830:    PETSc interface to MUMPS
1831: */
1832: #ifdef PETSC_HAVE_MUMPS
1833: PETSC_EXTERN PetscErrorCode MatMumpsSetIcntl(Mat,PetscInt,PetscInt);
1834: PETSC_EXTERN PetscErrorCode MatMumpsGetIcntl(Mat,PetscInt,PetscInt*);
1835: PETSC_EXTERN PetscErrorCode MatMumpsSetCntl(Mat,PetscInt,PetscReal);
1836: PETSC_EXTERN PetscErrorCode MatMumpsGetCntl(Mat,PetscInt,PetscReal*);

1838: PETSC_EXTERN PetscErrorCode MatMumpsGetInfo(Mat,PetscInt,PetscInt*);
1839: PETSC_EXTERN PetscErrorCode MatMumpsGetInfog(Mat,PetscInt,PetscInt*);
1840: PETSC_EXTERN PetscErrorCode MatMumpsGetRinfo(Mat,PetscInt,PetscReal*);
1841: PETSC_EXTERN PetscErrorCode MatMumpsGetRinfog(Mat,PetscInt,PetscReal*);
1842: PETSC_EXTERN PetscErrorCode MatMumpsGetInverse(Mat,Mat);
1843: PETSC_EXTERN PetscErrorCode MatMumpsGetInverseTranspose(Mat,Mat);
1844: #endif

1846: /*
1847:    PETSc interface to Mkl_Pardiso
1848: */
1849: #ifdef PETSC_HAVE_MKL_PARDISO
1850: PETSC_EXTERN PetscErrorCode MatMkl_PardisoSetCntl(Mat,PetscInt,PetscInt);
1851: #endif

1853: /*
1854:    PETSc interface to Mkl_CPardiso
1855: */
1856: #ifdef PETSC_HAVE_MKL_CPARDISO
1857: PETSC_EXTERN PetscErrorCode MatMkl_CPardisoSetCntl(Mat,PetscInt,PetscInt);
1858: #endif

1860: /*
1861:    PETSc interface to SUPERLU
1862: */
1863: #ifdef PETSC_HAVE_SUPERLU
1864: PETSC_EXTERN PetscErrorCode MatSuperluSetILUDropTol(Mat,PetscReal);
1865: #endif

1867: /*
1868:    PETSc interface to SUPERLU_DIST
1869: */
1870: #ifdef PETSC_HAVE_SUPERLU_DIST
1871: PETSC_EXTERN PetscErrorCode MatSuperluDistGetDiagU(Mat,PetscScalar*);
1872: #endif

1874: /*
1875:    PETSc interface to STRUMPACK
1876: */
1877: #ifdef PETSC_HAVE_STRUMPACK
1878: /*E
1879:     MatSTRUMPACKReordering - sparsity reducing ordering to be used in STRUMPACK

1881:     Level: intermediate
1882: E*/
1883: typedef enum {MAT_STRUMPACK_NATURAL=0,
1884:               MAT_STRUMPACK_METIS=1,
1885:               MAT_STRUMPACK_PARMETIS=2,
1886:               MAT_STRUMPACK_SCOTCH=3,
1887:               MAT_STRUMPACK_PTSCOTCH=4,
1888:               MAT_STRUMPACK_RCM=5} MatSTRUMPACKReordering;
1889: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetReordering(Mat,MatSTRUMPACKReordering);
1890: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetColPerm(Mat,PetscBool);
1891: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSRelTol(Mat,PetscReal);
1892: PETSC_DEPRECATED_FUNCTION("Use MatSTRUMPACKSetHSSRelTol() (since version 3.9)") PETSC_STATIC_INLINE PetscErrorCode MatSTRUMPACKSetHSSRelCompTol(Mat mat,PetscReal rtol) {return MatSTRUMPACKSetHSSRelTol(mat,rtol);}
1893: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSAbsTol(Mat,PetscReal);
1894: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSMinSepSize(Mat,PetscInt);
1895: PETSC_DEPRECATED_FUNCTION("Use MatSTRUMPACKSetHSSMinSepSize() (since version 3.9)") PETSC_STATIC_INLINE PetscErrorCode MatSTRUMPACKSetHSSMinSize(Mat mat,PetscInt hssminsize) {return MatSTRUMPACKSetHSSMinSepSize(mat,hssminsize);}
1896: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSMaxRank(Mat,PetscInt);
1897: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSLeafSize(Mat,PetscInt);
1898: #endif


1901: PETSC_EXTERN PetscErrorCode MatBindToCPU(Mat,PetscBool);
1902: PETSC_DEPRECATED_FUNCTION("Use MatBindToCPU (since v3.13)") PETSC_STATIC_INLINE PetscErrorCode MatPinToCPU(Mat A,PetscBool flg) {return MatBindToCPU(A,flg);}

1904: typedef struct _p_SplitCSRMat PetscSplitCSRDataStructure;

1906: #ifdef PETSC_HAVE_KOKKOS_KERNELS
1907: PETSC_EXTERN PetscErrorCode MatKokkosGetDeviceMatWrite(Mat,PetscSplitCSRDataStructure**);
1908: PETSC_EXTERN PetscErrorCode MatSeqAIJKokkosSetDeviceMat(Mat, PetscSplitCSRDataStructure *);
1909: PETSC_EXTERN PetscErrorCode MatSeqAIJKokkosGetDeviceMat(Mat, PetscSplitCSRDataStructure **);
1910: #endif

1912: #ifdef PETSC_HAVE_CUDA
1913: /*E
1914:     MatCUSPARSEStorageFormat - indicates the storage format for CUSPARSE (GPU)
1915:     matrices.

1917:     Not Collective

1919: +   MAT_CUSPARSE_CSR - Compressed Sparse Row
1920: .   MAT_CUSPARSE_ELL - Ellpack (requires CUDA 4.2 or later).
1921: -   MAT_CUSPARSE_HYB - Hybrid, a combination of Ellpack and Coordinate format (requires CUDA 4.2 or later).

1923:     Level: intermediate

1925:    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h

1927: .seealso: MatCUSPARSESetFormat(), MatCUSPARSEFormatOperation
1928: E*/

1930: typedef enum {MAT_CUSPARSE_CSR, MAT_CUSPARSE_ELL, MAT_CUSPARSE_HYB} MatCUSPARSEStorageFormat;

1932: /* these will be strings associated with enumerated type defined above */
1933: PETSC_EXTERN const char *const MatCUSPARSEStorageFormats[];

1935: /*E
1936:     MatCUSPARSEFormatOperation - indicates the operation of CUSPARSE (GPU)
1937:     matrices whose operation should use a particular storage format.

1939:     Not Collective

1941: +   MAT_CUSPARSE_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult
1942: .   MAT_CUSPARSE_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult
1943: .   MAT_CUSPARSE_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult
1944: -   MAT_CUSPARSE_ALL - sets the storage format for all CUSPARSE (GPU) matrices

1946:     Level: intermediate

1948: .seealso: MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat
1949: E*/
1950: typedef enum {MAT_CUSPARSE_MULT_DIAG, MAT_CUSPARSE_MULT_OFFDIAG, MAT_CUSPARSE_MULT, MAT_CUSPARSE_ALL} MatCUSPARSEFormatOperation;

1952: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
1953: PETSC_EXTERN PetscErrorCode MatCreateAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1954: PETSC_EXTERN PetscErrorCode MatCUSPARSESetFormat(Mat,MatCUSPARSEFormatOperation,MatCUSPARSEStorageFormat);
1955: typedef struct Mat_SeqAIJCUSPARSETriFactors* Mat_SeqAIJCUSPARSETriFactors_p;
1956: PETSC_INTERN PetscErrorCode MatSeqAIJCUSPARSETriFactors_Reset(Mat_SeqAIJCUSPARSETriFactors_p*);

1958: PETSC_EXTERN PetscErrorCode MatCUSPARSEGetDeviceMatWrite(Mat,PetscSplitCSRDataStructure**);

1960: PETSC_EXTERN PetscErrorCode MatCreateDenseCUDA(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
1961: PETSC_EXTERN PetscErrorCode MatCreateSeqDenseCUDA(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
1962: PETSC_EXTERN PetscErrorCode MatMPIDenseCUDASetPreallocation(Mat,PetscScalar[]);
1963: PETSC_EXTERN PetscErrorCode MatSeqDenseCUDASetPreallocation(Mat,PetscScalar[]);
1964: PETSC_EXTERN PetscErrorCode MatDenseCUDAGetArrayWrite(Mat,PetscScalar**);
1965: PETSC_EXTERN PetscErrorCode MatDenseCUDAGetArrayRead(Mat,const PetscScalar**);
1966: PETSC_EXTERN PetscErrorCode MatDenseCUDAGetArray(Mat,PetscScalar**);
1967: PETSC_EXTERN PetscErrorCode MatDenseCUDARestoreArrayWrite(Mat,PetscScalar**);
1968: PETSC_EXTERN PetscErrorCode MatDenseCUDARestoreArrayRead(Mat,const PetscScalar**);
1969: PETSC_EXTERN PetscErrorCode MatDenseCUDARestoreArray(Mat,PetscScalar**);
1970: PETSC_EXTERN PetscErrorCode MatDenseCUDAPlaceArray(Mat,const PetscScalar*);
1971: PETSC_EXTERN PetscErrorCode MatDenseCUDAReplaceArray(Mat,const PetscScalar*);
1972: PETSC_EXTERN PetscErrorCode MatDenseCUDAResetArray(Mat);

1974: #endif

1976: #if defined(PETSC_HAVE_VIENNACL)
1977: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
1978: PETSC_EXTERN PetscErrorCode MatCreateAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1979: #endif

1981: /*
1982:    PETSc interface to FFTW
1983: */
1984: #if defined(PETSC_HAVE_FFTW)
1985: PETSC_EXTERN PetscErrorCode VecScatterPetscToFFTW(Mat,Vec,Vec);
1986: PETSC_EXTERN PetscErrorCode VecScatterFFTWToPetsc(Mat,Vec,Vec);
1987: PETSC_EXTERN PetscErrorCode MatCreateVecsFFTW(Mat,Vec*,Vec*,Vec*);
1988: #endif

1990: #if defined(PETSC_HAVE_SCALAPACK)
1991: PETSC_EXTERN PetscErrorCode MatCreateScaLAPACK(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
1992: PETSC_EXTERN PetscErrorCode MatScaLAPACKSetBlockSizes(Mat,PetscInt,PetscInt);
1993: PETSC_EXTERN PetscErrorCode MatScaLAPACKGetBlockSizes(Mat,PetscInt*,PetscInt*);
1994: #endif

1996: PETSC_EXTERN PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*);
1997: PETSC_EXTERN PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*);
1998: PETSC_EXTERN PetscErrorCode MatNestGetISs(Mat,IS[],IS[]);
1999: PETSC_EXTERN PetscErrorCode MatNestGetLocalISs(Mat,IS[],IS[]);
2000: PETSC_EXTERN PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***);
2001: PETSC_EXTERN PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*);
2002: PETSC_EXTERN PetscErrorCode MatNestSetVecType(Mat,VecType);
2003: PETSC_EXTERN PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]);
2004: PETSC_EXTERN PetscErrorCode MatNestSetSubMat(Mat,PetscInt,PetscInt,Mat);

2006: PETSC_EXTERN PetscErrorCode MatChop(Mat,PetscReal);
2007: PETSC_EXTERN PetscErrorCode MatComputeBandwidth(Mat,PetscReal,PetscInt*);

2009: PETSC_EXTERN PetscErrorCode MatSubdomainsCreateCoalesce(Mat,PetscInt,PetscInt*,IS**);

2011: PETSC_EXTERN PetscErrorCode MatPreallocatorPreallocate(Mat,PetscBool,Mat);

2013: PETSC_INTERN PetscErrorCode MatHeaderMerge(Mat,Mat*);
2014: PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat,Mat*);

2016: #endif