Actual source code: petscmat.h

  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 MatFactorGetUseOrdering(Mat, PetscBool*);
164: PETSC_EXTERN PetscErrorCode MatFactorGetSolverType(Mat,MatSolverType*);
165: PETSC_EXTERN PetscErrorCode MatGetFactorType(Mat,MatFactorType*);
166: PETSC_EXTERN PetscErrorCode MatSetFactorType(Mat,MatFactorType);
167: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*MatSolverFunction)(Mat,MatFactorType,Mat*);
168: PETSC_EXTERN PetscErrorCode MatSolverTypeRegister(MatSolverType,MatType,MatFactorType,MatSolverFunction);
169: PETSC_EXTERN PetscErrorCode MatSolverTypeGet(MatSolverType,MatType,MatFactorType,PetscBool*,PetscBool*,MatSolverFunction*);
170: typedef MatSolverType MatSolverPackage PETSC_DEPRECATED_TYPEDEF("Use MatSolverType (since version 3.9)");
171: PETSC_DEPRECATED_FUNCTION("Use MatSolverTypeRegister() (since version 3.9)") PETSC_STATIC_INLINE PetscErrorCode MatSolverPackageRegister(MatSolverType stype,MatType mtype,MatFactorType ftype,MatSolverFunction f)
172: { return MatSolverTypeRegister(stype,mtype,ftype,f); }
173: 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)
174: { return MatSolverTypeGet(stype,mtype,ftype,foundmtype,foundstype,f); }

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

179:     Level: beginner

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

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

189:    Level: beginner

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

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

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

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

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

228:     Level: beginner

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

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

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

240:     Level: beginner

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

246: PETSC_EXTERN PetscErrorCode MatInitializePackage(void);

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

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

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

269:     Level: beginner

271: $  SAME_NONZERO_PATTERN  - the two matrices have identical nonzero patterns
272: $  DIFFERENT_NONZERO_PATTERN - the two matrices may have different nonzero patterns
273: $  SUBSET_NONZERO_PATTERN - the nonzero pattern of the second matrix is a subset of the nonzero pattern of the first matrix
274: $  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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

392:    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).
393:    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.

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

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

399:    Level: beginner

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

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

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

415:     Level: beginner

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



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

429:     Level: beginner

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

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

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

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

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

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

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

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

541:     Level: beginner

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

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

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

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

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


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

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

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

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

586:    Level: intermediate

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

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

604:     Level: beginner

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

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

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

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

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

655: PETSC_EXTERN PetscErrorCode MatCreateSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
656: 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);}
657: PETSC_EXTERN PetscErrorCode MatCreateSubMatricesMPI(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
658: 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);}
659: PETSC_EXTERN PetscErrorCode MatDestroyMatrices(PetscInt,Mat *[]);
660: PETSC_EXTERN PetscErrorCode MatDestroySubMatrices(PetscInt,Mat *[]);
661: PETSC_EXTERN PetscErrorCode MatCreateSubMatrix(Mat,IS,IS,MatReuse,Mat *);
662: 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);}
663: PETSC_EXTERN PetscErrorCode MatGetLocalSubMatrix(Mat,IS,IS,Mat*);
664: PETSC_EXTERN PetscErrorCode MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*);
665: PETSC_EXTERN PetscErrorCode MatGetSeqNonzeroStructure(Mat,Mat*);
666: PETSC_EXTERN PetscErrorCode MatDestroySeqNonzeroStructure(Mat*);

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

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

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

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

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

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

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

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

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

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

713: PETSC_EXTERN PetscErrorCode MatInterpolate(Mat,Vec,Vec);
714: PETSC_EXTERN PetscErrorCode MatInterpolateAdd(Mat,Vec,Vec,Vec);
715: PETSC_EXTERN PetscErrorCode MatRestrict(Mat,Vec,Vec);
716: PETSC_EXTERN PetscErrorCode MatMatInterpolate(Mat,Mat,Mat*);
717: PETSC_EXTERN PetscErrorCode MatMatInterpolateAdd(Mat,Mat,Mat,Mat*);
718: PETSC_EXTERN PetscErrorCode MatMatRestrict(Mat,Mat,Mat*);
719: PETSC_EXTERN PetscErrorCode MatCreateVecs(Mat,Vec*,Vec*);
720: 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);}
721: PETSC_EXTERN PetscErrorCode MatCreateRedundantMatrix(Mat,PetscInt,MPI_Comm,MatReuse,Mat*);
722: PETSC_EXTERN PetscErrorCode MatGetMultiProcBlock(Mat,MPI_Comm,MatReuse,Mat*);
723: PETSC_EXTERN PetscErrorCode MatFindZeroDiagonals(Mat,IS*);
724: PETSC_EXTERN PetscErrorCode MatFindOffBlockDiagonalEntries(Mat,IS*);
725: PETSC_EXTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);

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

730:    Not collective

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

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

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

747:    Level: beginner

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

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

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

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

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

765:    Collective

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

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

776:    Level: intermediate

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

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

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

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

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

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

805:    Not Collective

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

817:    Level: intermediate

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

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

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

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

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

845:    Not Collective

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

857:    Level: intermediate

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

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

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

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

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

886:    Not Collective

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

898:    Level: intermediate

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

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

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

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

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

926:    Not Collective

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

937:    Level: intermediate

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

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

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

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

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

965:    Not Collective

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

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

976:    Level: intermediate

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

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

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

985: .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock()
986:           MatPreallocateInitialize(), MatPreallocateSetLocal()
987: M*/
988: #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
989: do { PetscInt __i; \
990:   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);\
991:   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);\
992:   for (__i=0; __i<nc; __i++) {\
993:     if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \
994:     else if (dnz[row - __rstart] < __ncols) dnz[row - __rstart]++;\
995:   }\
996: } while (0)

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

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

1006:    Not Collective

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

1016:    Level: intermediate

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

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

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

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

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

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

1043:    Not Collective

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

1053:    Level: intermediate

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

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

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

1062: .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
1063:           MatPreallocateSymmetricSetLocalBlock()
1064: M*/
1065: #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)


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

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

1076:    Collective

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

1082:    Level: intermediate

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

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

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

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

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

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

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

1109: #define MAT_SKIP_ALLOCATION -4

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

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

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

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

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

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

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

1151:    Level: beginner

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

1156: .seealso: MatGetOrdering()
1157: J*/
1158: typedef const char* MatOrderingType;
1159: #define MATORDERINGNATURAL        "natural"
1160: #define MATORDERINGND             "nd"
1161: #define MATORDERING1WD            "1wd"
1162: #define MATORDERINGRCM            "rcm"
1163: #define MATORDERINGQMD            "qmd"
1164: #define MATORDERINGROWLENGTH      "rowlength"
1165: #define MATORDERINGWBM            "wbm"
1166: #define MATORDERINGSPECTRAL       "spectral"
1167: #define MATORDERINGAMD            "amd"            /* only works if UMFPACK is installed with PETSc */
1168: #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 */
1169: #define MATORDERINGEXTERNAL       "external"       /* uses an ordering type internal to the factorization package */

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

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

1179: /*S
1180:     MatFactorShiftType - Numeric Shift for factorizations

1182:    Level: beginner

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

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

1193:     Level: beginner

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

1198: .seealso: MatGetFactor()
1199: S*/
1200: typedef enum {MAT_FACTOR_NOERROR,MAT_FACTOR_STRUCT_ZEROPIVOT,MAT_FACTOR_NUMERIC_ZEROPIVOT,MAT_FACTOR_OUTMEMORY,MAT_FACTOR_OTHER} MatFactorError;

1202: PETSC_EXTERN PetscErrorCode MatFactorGetError(Mat,MatFactorError*);
1203: PETSC_EXTERN PetscErrorCode MatFactorClearError(Mat);
1204: PETSC_EXTERN PetscErrorCode MatFactorGetErrorZeroPivot(Mat,PetscReal*,PetscInt*);

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

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

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

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

1217:    Level: developer

1219: .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
1220:           MatFactorInfoInitialize()

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

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

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

1272: /*E
1273:     MatSORType - What type of (S)SOR to perform

1275:     Level: beginner

1277:    May be bitwise ORd together

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

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

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

1291: /*
1292:     These routines are for efficiently computing Jacobians via finite differences.
1293: */

1295: /*S
1296:      MatColoring - Object for managing the coloring of matrices.

1298:    Level: beginner

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

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

1307: .seealso:  MatFDColoringCreate(), MatColoringWeightType, ISColoring, MatFDColoring, DMCreateColoring(), MatColoringCreate(), MatOrdering, MatPartitioning
1308: S*/
1309: typedef struct _p_MatColoring* MatColoring;

1311: /*J
1312:     MatColoringType - String with the name of a PETSc matrix coloring

1314:    Level: beginner

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

1327: /*E
1328:    MatColoringWeightType - Type of weight scheme

1330:     Not Collective

1332: +   MAT_COLORING_RANDOM  - Random weights
1333: .   MAT_COLORING_LEXICAL - Lexical weighting based upon global numbering.
1334: -   MAT_COLORING_LF      - Last-first weighting.

1336:     Level: intermediate

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

1340: .seealso: MatColoring, MatColoringCreate()
1341: E*/
1342: typedef enum {MAT_COLORING_WEIGHT_RANDOM,MAT_COLORING_WEIGHT_LEXICAL,MAT_COLORING_WEIGHT_LF,MAT_COLORING_WEIGHT_SL} MatColoringWeightType;

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

1364: /*S
1365:      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1366:         and coloring

1368:    Level: beginner

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

1373: .seealso:  MatFDColoringCreate(), MatColoring, DMCreateColoring()
1374: S*/
1375: typedef struct _p_MatFDColoring* MatFDColoring;

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

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

1394:    Level: beginner

1396: .seealso:  MatTransposeColoringCreate()
1397: S*/
1398: typedef struct _p_MatTransposeColoring* MatTransposeColoring;

1400: PETSC_EXTERN PetscErrorCode MatTransposeColoringCreate(Mat,ISColoring,MatTransposeColoring *);
1401: PETSC_EXTERN PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring,Mat,Mat);
1402: PETSC_EXTERN PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring,Mat,Mat);
1403: PETSC_EXTERN PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring*);

1405: /*
1406:     These routines are for partitioning matrices: currently used only
1407:   for adjacency matrix, MatCreateMPIAdj().
1408: */

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

1413:    Level: beginner

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

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

1422: .seealso:  MatPartitioningCreate(), MatPartitioningType, MatColoring, MatGetOrdering()
1423: S*/
1424: typedef struct _p_MatPartitioning* MatPartitioning;

1426: /*J
1427:     MatPartitioningType - String with the name of a PETSc matrix partitioning

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

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

1462: PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetRepartition(MatPartitioning part);
1463: PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
1464: PETSC_EXTERN PetscErrorCode MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *);

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

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

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

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

1505: PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning,PetscReal);
1506: PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning,PetscReal*);
1507: PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning,MPPTScotchStrategyType);
1508: PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning,MPPTScotchStrategyType*);

1510: /*
1511:  * hierarchical partitioning
1512:  */
1513: PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetFineparts(MatPartitioning,IS*);
1514: PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetCoarseparts(MatPartitioning,IS*);
1515: PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNcoarseparts(MatPartitioning,PetscInt);
1516: PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNfineparts(MatPartitioning, PetscInt);

1518: PETSC_EXTERN PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*);
1519: PETSC_EXTERN PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*);

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

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

1698: PETSC_EXTERN PetscErrorCode MatMPIBAIJSetHashTableFactor(Mat,PetscReal);

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

1709: /*S
1710:      MatNullSpace - Object that removes a null space from a vector, i.e.
1711:          orthogonalizes the vector to a subspace

1713:    Level: advanced

1715:   Users manual sections:
1716: .   sec_singular

1718: .seealso:  MatNullSpaceCreate()
1719: S*/
1720: typedef struct _p_MatNullSpace* MatNullSpace;

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

1737: PETSC_EXTERN PetscErrorCode MatReorderingSeqSBAIJ(Mat,IS);
1738: PETSC_EXTERN PetscErrorCode MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
1739: PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);

1741: PETSC_EXTERN PetscErrorCode MatCreateMAIJ(Mat,PetscInt,Mat*);
1742: PETSC_EXTERN PetscErrorCode MatMAIJRedimension(Mat,PetscInt,Mat*);
1743: PETSC_EXTERN PetscErrorCode MatMAIJGetAIJ(Mat,Mat*);

1745: PETSC_EXTERN PetscErrorCode MatComputeOperator(Mat,MatType,Mat*);
1746: PETSC_EXTERN PetscErrorCode MatComputeOperatorTranspose(Mat,MatType,Mat*);

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

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

1766: PETSC_EXTERN PetscErrorCode MatDiagonalScaleLocal(Mat,Vec);

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

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

1786:     Notes:
1787:     MATMFFD is a specific MatType which uses the MatMFFD data structure

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

1791:     Level: developer

1793: .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister()
1794: S*/
1795: typedef struct _p_MatMFFD* MatMFFD;

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

1800:    Level: beginner

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

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

1811: PETSC_EXTERN PetscErrorCode MatMFFDDSSetUmin(Mat,PetscReal);
1812: PETSC_EXTERN PetscErrorCode MatMFFDWPSetComputeNormU(Mat,PetscBool);

1814: PETSC_EXTERN PetscErrorCode MatFDColoringSetType(MatFDColoring,MatMFFDType);

1816: PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
1817: PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);

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

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

1835: PETSC_EXTERN PetscErrorCode MatMumpsGetInfo(Mat,PetscInt,PetscInt*);
1836: PETSC_EXTERN PetscErrorCode MatMumpsGetInfog(Mat,PetscInt,PetscInt*);
1837: PETSC_EXTERN PetscErrorCode MatMumpsGetRinfo(Mat,PetscInt,PetscReal*);
1838: PETSC_EXTERN PetscErrorCode MatMumpsGetRinfog(Mat,PetscInt,PetscReal*);
1839: PETSC_EXTERN PetscErrorCode MatMumpsGetInverse(Mat,Mat);
1840: PETSC_EXTERN PetscErrorCode MatMumpsGetInverseTranspose(Mat,Mat);
1841: #endif

1843: /*
1844:    PETSc interface to Mkl_Pardiso
1845: */
1846: #ifdef PETSC_HAVE_MKL_PARDISO
1847: PETSC_EXTERN PetscErrorCode MatMkl_PardisoSetCntl(Mat,PetscInt,PetscInt);
1848: #endif

1850: /*
1851:    PETSc interface to Mkl_CPardiso
1852: */
1853: #ifdef PETSC_HAVE_MKL_CPARDISO
1854: PETSC_EXTERN PetscErrorCode MatMkl_CPardisoSetCntl(Mat,PetscInt,PetscInt);
1855: #endif

1857: /*
1858:    PETSc interface to SUPERLU
1859: */
1860: #ifdef PETSC_HAVE_SUPERLU
1861: PETSC_EXTERN PetscErrorCode MatSuperluSetILUDropTol(Mat,PetscReal);
1862: #endif

1864: /*
1865:    PETSc interface to SUPERLU_DIST
1866: */
1867: #ifdef PETSC_HAVE_SUPERLU_DIST
1868: PETSC_EXTERN PetscErrorCode MatSuperluDistGetDiagU(Mat,PetscScalar*);
1869: #endif

1871: /*
1872:    PETSc interface to STRUMPACK
1873: */
1874: #ifdef PETSC_HAVE_STRUMPACK
1875: /*E
1876:     MatSTRUMPACKReordering - sparsity reducing ordering to be used in STRUMPACK

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


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

1901: typedef struct _p_SplitCSRMat PetscSplitCSRDataStructure;

1903: #ifdef PETSC_HAVE_KOKKOS_KERNELS
1904: PETSC_EXTERN PetscErrorCode MatKokkosGetDeviceMatWrite(Mat,PetscSplitCSRDataStructure**);
1905: PETSC_EXTERN PetscErrorCode MatSeqAIJKokkosSetDeviceMat(Mat, PetscSplitCSRDataStructure *);
1906: PETSC_EXTERN PetscErrorCode MatSeqAIJKokkosGetDeviceMat(Mat, PetscSplitCSRDataStructure **);
1907: #endif

1909: #ifdef PETSC_HAVE_CUDA
1910: /*E
1911:     MatCUSPARSEStorageFormat - indicates the storage format for CUSPARSE (GPU)
1912:     matrices.

1914:     Not Collective

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

1920:     Level: intermediate

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

1924: .seealso: MatCUSPARSESetFormat(), MatCUSPARSEFormatOperation
1925: E*/

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

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

1932: /*E
1933:     MatCUSPARSEFormatOperation - indicates the operation of CUSPARSE (GPU)
1934:     matrices whose operation should use a particular storage format.

1936:     Not Collective

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

1943:     Level: intermediate

1945: .seealso: MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat
1946: E*/
1947: typedef enum {MAT_CUSPARSE_MULT_DIAG, MAT_CUSPARSE_MULT_OFFDIAG, MAT_CUSPARSE_MULT, MAT_CUSPARSE_ALL} MatCUSPARSEFormatOperation;

1949: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
1950: PETSC_EXTERN PetscErrorCode MatCreateAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1951: PETSC_EXTERN PetscErrorCode MatCUSPARSESetFormat(Mat,MatCUSPARSEFormatOperation,MatCUSPARSEStorageFormat);

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

1955: PETSC_EXTERN PetscErrorCode MatCreateDenseCUDA(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
1956: PETSC_EXTERN PetscErrorCode MatCreateSeqDenseCUDA(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
1957: PETSC_EXTERN PetscErrorCode MatMPIDenseCUDASetPreallocation(Mat,PetscScalar[]);
1958: PETSC_EXTERN PetscErrorCode MatSeqDenseCUDASetPreallocation(Mat,PetscScalar[]);
1959: PETSC_EXTERN PetscErrorCode MatDenseCUDAGetArrayWrite(Mat,PetscScalar**);
1960: PETSC_EXTERN PetscErrorCode MatDenseCUDAGetArrayRead(Mat,const PetscScalar**);
1961: PETSC_EXTERN PetscErrorCode MatDenseCUDAGetArray(Mat,PetscScalar**);
1962: PETSC_EXTERN PetscErrorCode MatDenseCUDARestoreArrayWrite(Mat,PetscScalar**);
1963: PETSC_EXTERN PetscErrorCode MatDenseCUDARestoreArrayRead(Mat,const PetscScalar**);
1964: PETSC_EXTERN PetscErrorCode MatDenseCUDARestoreArray(Mat,PetscScalar**);
1965: PETSC_EXTERN PetscErrorCode MatDenseCUDAPlaceArray(Mat,const PetscScalar*);
1966: PETSC_EXTERN PetscErrorCode MatDenseCUDAReplaceArray(Mat,const PetscScalar*);
1967: PETSC_EXTERN PetscErrorCode MatDenseCUDAResetArray(Mat);

1969: #endif

1971: #if defined(PETSC_HAVE_VIENNACL)
1972: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
1973: PETSC_EXTERN PetscErrorCode MatCreateAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1974: #endif

1976: /*
1977:    PETSc interface to FFTW
1978: */
1979: #if defined(PETSC_HAVE_FFTW)
1980: PETSC_EXTERN PetscErrorCode VecScatterPetscToFFTW(Mat,Vec,Vec);
1981: PETSC_EXTERN PetscErrorCode VecScatterFFTWToPetsc(Mat,Vec,Vec);
1982: PETSC_EXTERN PetscErrorCode MatCreateVecsFFTW(Mat,Vec*,Vec*,Vec*);
1983: #endif

1985: #if defined(PETSC_HAVE_SCALAPACK)
1986: PETSC_EXTERN PetscErrorCode MatCreateScaLAPACK(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
1987: PETSC_EXTERN PetscErrorCode MatScaLAPACKSetBlockSizes(Mat,PetscInt,PetscInt);
1988: PETSC_EXTERN PetscErrorCode MatScaLAPACKGetBlockSizes(Mat,PetscInt*,PetscInt*);
1989: #endif

1991: PETSC_EXTERN PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*);
1992: PETSC_EXTERN PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*);
1993: PETSC_EXTERN PetscErrorCode MatNestGetISs(Mat,IS[],IS[]);
1994: PETSC_EXTERN PetscErrorCode MatNestGetLocalISs(Mat,IS[],IS[]);
1995: PETSC_EXTERN PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***);
1996: PETSC_EXTERN PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*);
1997: PETSC_EXTERN PetscErrorCode MatNestSetVecType(Mat,VecType);
1998: PETSC_EXTERN PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]);
1999: PETSC_EXTERN PetscErrorCode MatNestSetSubMat(Mat,PetscInt,PetscInt,Mat);

2001: PETSC_EXTERN PetscErrorCode MatChop(Mat,PetscReal);
2002: PETSC_EXTERN PetscErrorCode MatComputeBandwidth(Mat,PetscReal,PetscInt*);

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

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

2008: PETSC_INTERN PetscErrorCode MatHeaderMerge(Mat,Mat*);
2009: PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat,Mat*);

2011: #endif