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 MATCENTERING       "centering"
 63: #define MATDENSE           "dense"
 64: #define MATDENSECUDA       "densecuda"
 65: #define MATSEQDENSE        "seqdense"
 66: #define MATSEQDENSECUDA    "seqdensecuda"
 67: #define MATMPIDENSE        "mpidense"
 68: #define MATMPIDENSECUDA    "mpidensecuda"
 69: #define MATELEMENTAL       "elemental"
 70: #define MATSCALAPACK       "scalapack"
 71: #define MATBAIJ            "baij"
 72: #define MATSEQBAIJ         "seqbaij"
 73: #define MATMPIBAIJ         "mpibaij"
 74: #define MATMPIADJ          "mpiadj"
 75: #define MATSBAIJ           "sbaij"
 76: #define MATSEQSBAIJ        "seqsbaij"
 77: #define MATMPISBAIJ        "mpisbaij"
 78: #define MATMFFD            "mffd"
 79: #define MATNORMAL          "normal"
 80: #define MATNORMALHERMITIAN "normalh"
 81: #define MATLRC             "lrc"
 82: #define MATSCATTER         "scatter"
 83: #define MATBLOCKMAT        "blockmat"
 84: #define MATCOMPOSITE       "composite"
 85: #define MATFFT             "fft"
 86: #define MATFFTW            "fftw"
 87: #define MATSEQCUFFT        "seqcufft"
 88: #define MATTRANSPOSEMAT    "transpose"
 89: #define MATSCHURCOMPLEMENT "schurcomplement"
 90: #define MATPYTHON          "python"
 91: #define MATHYPRE           "hypre"
 92: #define MATHYPRESTRUCT     "hyprestruct"
 93: #define MATHYPRESSTRUCT    "hypresstruct"
 94: #define MATSUBMATRIX       "submatrix"
 95: #define MATLOCALREF        "localref"
 96: #define MATNEST            "nest"
 97: #define MATPREALLOCATOR    "preallocator"
 98: #define MATSELL            "sell"
 99: #define MATSEQSELL         "seqsell"
100: #define MATMPISELL         "mpisell"
101: #define MATDUMMY           "dummy"
102: #define MATLMVM            "lmvm"
103: #define MATLMVMDFP         "lmvmdfp"
104: #define MATLMVMBFGS        "lmvmbfgs"
105: #define MATLMVMSR1         "lmvmsr1"
106: #define MATLMVMBROYDEN     "lmvmbroyden"
107: #define MATLMVMBADBROYDEN  "lmvmbadbroyden"
108: #define MATLMVMSYMBROYDEN  "lmvmsymbroyden"
109: #define MATLMVMSYMBADBROYDEN "lmvmsymbadbroyden"
110: #define MATLMVMDIAGBROYDEN   "lmvmdiagbroyden"
111: #define MATCONSTANTDIAGONAL  "constantdiagonal"
112: #define MATHTOOL           "htool"
113: #define MATH2OPUS          "h2opus"

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

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

120:    Level: beginner

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

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

152: /*E
153:     MatFactorType - indicates what type of factorization is requested

155:     Level: beginner

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

159: .seealso: MatSolverType, MatGetFactor(), MatGetFactorAvailable(), MatSolverTypeRegister()
160: E*/
161: 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;
162: PETSC_EXTERN const char *const MatFactorTypes[];

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

180: /*E
181:     MatProductType - indicates what type of matrix product is requested

183:     Level: beginner

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

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

193:    Level: beginner

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

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

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

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

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

232:     Level: beginner

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

236: .seealso: MatCreateSubMatrices(), MatCreateSubMatrix(), MatDestroyMatrices(), MatConvert()
237: E*/
238: typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX,MAT_IGNORE_MATRIX,MAT_INPLACE_MATRIX} MatReuse;

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

244:     Level: beginner

246: .seealso: MatGetSeqNonzeroStructure()
247: E*/
248: typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatCreateSubMatrixOption;

250: PETSC_EXTERN PetscErrorCode MatInitializePackage(void);

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

266: PETSC_EXTERN PetscFunctionList MatList;
267: PETSC_EXTERN PetscFunctionList MatColoringList;
268: PETSC_EXTERN PetscFunctionList MatPartitioningList;

270: /*E
271:     MatStructure - Indicates if two matrices have the same nonzero structure

273:     Level: beginner

275: $  SAME_NONZERO_PATTERN  - the two matrices have identical nonzero patterns
276: $  DIFFERENT_NONZERO_PATTERN - the two matrices may have different nonzero patterns
277: $  SUBSET_NONZERO_PATTERN - the nonzero pattern of the second matrix is a subset of the nonzero pattern of the first matrix
278: $  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

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

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

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

295: PETSC_EXTERN PetscErrorCode MatCreateSeqSELL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
296: PETSC_EXTERN PetscErrorCode MatCreateSELL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
297: PETSC_EXTERN PetscErrorCode MatSeqSELLSetPreallocation(Mat,PetscInt,const PetscInt[]);
298: PETSC_EXTERN PetscErrorCode MatMPISELLSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);

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

309: PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
310: PETSC_EXTERN PetscErrorCode MatCreateBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
311: PETSC_EXTERN PetscErrorCode MatCreateMPIBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat*);

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

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

319: PETSC_EXTERN PetscErrorCode MatCreateSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
320: PETSC_EXTERN PetscErrorCode MatCreateMPISBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
321: PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
322: PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
323: PETSC_EXTERN PetscErrorCode MatXAIJSetPreallocation(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscInt[],const PetscInt[]);

325: PETSC_EXTERN PetscErrorCode MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*);
326: PETSC_EXTERN PetscErrorCode MatCreateCentering(MPI_Comm,PetscInt,PetscInt,Mat*);
327: PETSC_EXTERN PetscErrorCode MatCreateNormal(Mat,Mat*);
328: PETSC_EXTERN PetscErrorCode MatCreateNormalHermitian(Mat,Mat*);
329: PETSC_EXTERN PetscErrorCode MatCreateLRC(Mat,Mat,Vec,Mat,Mat*);
330: PETSC_EXTERN PetscErrorCode MatLRCGetMats(Mat,Mat*,Mat*,Vec*,Mat*);
331: PETSC_EXTERN PetscErrorCode MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,ISLocalToGlobalMapping,Mat*);
332: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
333: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);

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

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

356: PETSC_EXTERN PetscErrorCode MatCreateTranspose(Mat,Mat*);
357: PETSC_EXTERN PetscErrorCode MatTransposeGetMat(Mat,Mat*);
358: PETSC_EXTERN PetscErrorCode MatCreateHermitianTranspose(Mat,Mat*);
359: PETSC_EXTERN PetscErrorCode MatHermitianTransposeGetMat(Mat,Mat*);
360: PETSC_EXTERN PetscErrorCode MatNormalGetMat(Mat,Mat*);
361: PETSC_EXTERN PetscErrorCode MatNormalHermitianGetMat(Mat,Mat*);
362: PETSC_EXTERN PetscErrorCode MatCreateSubMatrixVirtual(Mat,IS,IS,Mat*);
363: PETSC_EXTERN PetscErrorCode MatSubMatrixVirtualUpdate(Mat,Mat,IS,IS);
364: PETSC_EXTERN PetscErrorCode MatCreateLocalRef(Mat,IS,IS,Mat*);
365: PETSC_EXTERN PetscErrorCode MatCreateConstantDiagonal(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar,Mat*);

367: #if defined(PETSC_HAVE_HYPRE)
368: PETSC_EXTERN PetscErrorCode MatHYPRESetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
369: #endif

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

373: PETSC_EXTERN PetscErrorCode MatResetPreallocation(Mat);
374: PETSC_EXTERN PetscErrorCode MatSetUp(Mat);
375: PETSC_EXTERN PetscErrorCode MatDestroy(Mat*);
376: PETSC_EXTERN PetscErrorCode MatGetNonzeroState(Mat,PetscObjectState*);

378: PETSC_EXTERN PetscErrorCode MatConjugate(Mat);
379: PETSC_EXTERN PetscErrorCode MatRealPart(Mat);
380: PETSC_EXTERN PetscErrorCode MatImaginaryPart(Mat);
381: PETSC_EXTERN PetscErrorCode MatGetDiagonalBlock(Mat,Mat*);
382: PETSC_EXTERN PetscErrorCode MatGetTrace(Mat,PetscScalar*);
383: PETSC_EXTERN PetscErrorCode MatInvertBlockDiagonal(Mat,const PetscScalar **);
384: PETSC_EXTERN PetscErrorCode MatInvertVariableBlockDiagonal(Mat,PetscInt,const PetscInt*,PetscScalar*);
385: PETSC_EXTERN PetscErrorCode MatInvertBlockDiagonalMat(Mat,Mat);

387: /* ------------------------------------------------------------*/
388: PETSC_EXTERN PetscErrorCode MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
389: PETSC_EXTERN PetscErrorCode MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
390: PETSC_EXTERN PetscErrorCode MatSetValuesRow(Mat,PetscInt,const PetscScalar[]);
391: PETSC_EXTERN PetscErrorCode MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]);
392: PETSC_EXTERN PetscErrorCode MatSetValuesBatch(Mat,PetscInt,PetscInt,PetscInt[],const PetscScalar[]);
393: PETSC_EXTERN PetscErrorCode MatSetRandom(Mat,PetscRandom);

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

399:    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).
400:    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.

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

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

406:    Level: beginner

408: .seealso:  MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockedStencil(), DMDAVecGetArray(), DMDAVecGetArrayF90()
409: S*/
410: typedef struct {
411:   PetscInt k,j,i,c;
412: } MatStencil;

414: PETSC_EXTERN PetscErrorCode MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
415: PETSC_EXTERN PetscErrorCode MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
416: PETSC_EXTERN PetscErrorCode MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt);

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

422:     Level: beginner

424: .seealso: MatAssemblyBegin(), MatAssemblyEnd()
425: E*/
426: typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType;
427: PETSC_EXTERN PetscErrorCode MatAssemblyBegin(Mat,MatAssemblyType);
428: PETSC_EXTERN PetscErrorCode MatAssemblyEnd(Mat,MatAssemblyType);
429: PETSC_EXTERN PetscErrorCode MatAssembled(Mat,PetscBool *);

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

434:     Level: beginner

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

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

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

473: PETSC_EXTERN const char *const *MatOptions;
474: PETSC_EXTERN PetscErrorCode MatSetOption(Mat,MatOption,PetscBool);
475: PETSC_EXTERN PetscErrorCode MatGetOption(Mat,MatOption,PetscBool*);
476: PETSC_EXTERN PetscErrorCode MatPropagateSymmetryOptions(Mat,Mat);
477: PETSC_EXTERN PetscErrorCode MatGetType(Mat,MatType*);

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

516: PETSC_EXTERN PetscErrorCode MatDenseGetColumn(Mat,PetscInt,PetscScalar*[]);
517: PETSC_EXTERN PetscErrorCode MatDenseRestoreColumn(Mat,PetscScalar*[]);
518: PETSC_EXTERN PetscErrorCode MatDenseGetColumnVec(Mat,PetscInt,Vec*);
519: PETSC_EXTERN PetscErrorCode MatDenseRestoreColumnVec(Mat,PetscInt,Vec*);
520: PETSC_EXTERN PetscErrorCode MatDenseGetColumnVecRead(Mat,PetscInt,Vec*);
521: PETSC_EXTERN PetscErrorCode MatDenseRestoreColumnVecRead(Mat,PetscInt,Vec*);
522: PETSC_EXTERN PetscErrorCode MatDenseGetColumnVecWrite(Mat,PetscInt,Vec*);
523: PETSC_EXTERN PetscErrorCode MatDenseRestoreColumnVecWrite(Mat,PetscInt,Vec*);
524: PETSC_EXTERN PetscErrorCode MatDenseGetSubMatrix(Mat,PetscInt,PetscInt,Mat*);
525: PETSC_EXTERN PetscErrorCode MatDenseRestoreSubMatrix(Mat,Mat*);

527: PETSC_EXTERN PetscErrorCode MatMult(Mat,Vec,Vec);
528: PETSC_EXTERN PetscErrorCode MatMultDiagonalBlock(Mat,Vec,Vec);
529: PETSC_EXTERN PetscErrorCode MatMultAdd(Mat,Vec,Vec,Vec);
530: PETSC_EXTERN PetscErrorCode MatMultTranspose(Mat,Vec,Vec);
531: PETSC_EXTERN PetscErrorCode MatMultHermitianTranspose(Mat,Vec,Vec);
532: PETSC_EXTERN PetscErrorCode MatIsTranspose(Mat,Mat,PetscReal,PetscBool *);
533: PETSC_EXTERN PetscErrorCode MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *);
534: PETSC_EXTERN PetscErrorCode MatMultTransposeAdd(Mat,Vec,Vec,Vec);
535: PETSC_EXTERN PetscErrorCode MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec);
536: PETSC_EXTERN PetscErrorCode MatMultConstrained(Mat,Vec,Vec);
537: PETSC_EXTERN PetscErrorCode MatMultTransposeConstrained(Mat,Vec,Vec);
538: PETSC_EXTERN PetscErrorCode MatMatSolve(Mat,Mat,Mat);
539: PETSC_EXTERN PetscErrorCode MatMatSolveTranspose(Mat,Mat,Mat);
540: PETSC_EXTERN PetscErrorCode MatMatTransposeSolve(Mat,Mat,Mat);
541: PETSC_EXTERN PetscErrorCode MatResidual(Mat,Vec,Vec,Vec);

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

547:     Level: beginner

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

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

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

564: .seealso: MatDuplicate()
565: E*/
566: typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption;

568: PETSC_EXTERN PetscErrorCode MatConvert(Mat,MatType,MatReuse,Mat*);
569: PETSC_EXTERN PetscErrorCode MatDuplicate(Mat,MatDuplicateOption,Mat*);

571: PETSC_EXTERN PetscErrorCode MatCopy(Mat,Mat,MatStructure);
572: PETSC_EXTERN PetscErrorCode MatView(Mat,PetscViewer);
573: PETSC_EXTERN PetscErrorCode MatIsSymmetric(Mat,PetscReal,PetscBool *);
574: PETSC_EXTERN PetscErrorCode MatIsStructurallySymmetric(Mat,PetscBool *);
575: PETSC_EXTERN PetscErrorCode MatIsHermitian(Mat,PetscReal,PetscBool *);
576: PETSC_EXTERN PetscErrorCode MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *);
577: PETSC_EXTERN PetscErrorCode MatIsHermitianKnown(Mat,PetscBool *,PetscBool *);
578: PETSC_EXTERN PetscErrorCode MatMissingDiagonal(Mat,PetscBool  *,PetscInt *);
579: PETSC_EXTERN PetscErrorCode MatLoad(Mat, PetscViewer);

581: PETSC_EXTERN PetscErrorCode MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
582: PETSC_EXTERN PetscErrorCode MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool  *);
583: PETSC_EXTERN PetscErrorCode MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
584: PETSC_EXTERN PetscErrorCode MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool  *);

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

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

591:    Level: intermediate

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

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

609:     Level: beginner

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

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

629: PETSC_EXTERN PetscErrorCode MatEqual(Mat,Mat,PetscBool*);
630: PETSC_EXTERN PetscErrorCode MatMultEqual(Mat,Mat,PetscInt,PetscBool*);
631: PETSC_EXTERN PetscErrorCode MatMultAddEqual(Mat,Mat,PetscInt,PetscBool*);
632: PETSC_EXTERN PetscErrorCode MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool*);
633: PETSC_EXTERN PetscErrorCode MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool*);
634: PETSC_EXTERN PetscErrorCode MatMatMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*);
635: PETSC_EXTERN PetscErrorCode MatTransposeMatMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*);
636: PETSC_EXTERN PetscErrorCode MatMatTransposeMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*);
637: PETSC_EXTERN PetscErrorCode MatPtAPMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*);
638: PETSC_EXTERN PetscErrorCode MatRARtMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*);
639: PETSC_EXTERN PetscErrorCode MatIsLinear(Mat,PetscInt,PetscBool*);

641: PETSC_EXTERN PetscErrorCode MatNorm(Mat,NormType,PetscReal*);
642: PETSC_EXTERN PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal*);
643: PETSC_EXTERN PetscErrorCode MatGetColumnSums(Mat,PetscScalar*);
644: PETSC_EXTERN PetscErrorCode MatGetColumnSumsRealPart(Mat,PetscReal*);
645: PETSC_EXTERN PetscErrorCode MatGetColumnSumsImaginaryPart(Mat,PetscReal*);
646: PETSC_EXTERN PetscErrorCode MatGetColumnMeans(Mat,PetscScalar*);
647: PETSC_EXTERN PetscErrorCode MatGetColumnMeansRealPart(Mat,PetscReal*);
648: PETSC_EXTERN PetscErrorCode MatGetColumnMeansImaginaryPart(Mat,PetscReal*);
649: PETSC_EXTERN PetscErrorCode MatGetColumnReductions(Mat,PetscInt,PetscReal*);
650: PETSC_EXTERN PetscErrorCode MatZeroEntries(Mat);
651: PETSC_EXTERN PetscErrorCode MatSetInf(Mat);
652: PETSC_EXTERN PetscErrorCode MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
653: PETSC_EXTERN PetscErrorCode MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec);
654: PETSC_EXTERN PetscErrorCode MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec);
655: PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsStencil(Mat,PetscInt,const MatStencil[],PetscScalar,Vec,Vec);
656: PETSC_EXTERN PetscErrorCode MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
657: PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec);

659: PETSC_EXTERN PetscErrorCode MatGetSize(Mat,PetscInt*,PetscInt*);
660: PETSC_EXTERN PetscErrorCode MatGetLocalSize(Mat,PetscInt*,PetscInt*);
661: PETSC_EXTERN PetscErrorCode MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
662: PETSC_EXTERN PetscErrorCode MatGetOwnershipRanges(Mat,const PetscInt**);
663: PETSC_EXTERN PetscErrorCode MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*);
664: PETSC_EXTERN PetscErrorCode MatGetOwnershipRangesColumn(Mat,const PetscInt**);
665: PETSC_EXTERN PetscErrorCode MatGetOwnershipIS(Mat,IS*,IS*);

667: PETSC_EXTERN PetscErrorCode MatCreateSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
668: 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);}
669: PETSC_EXTERN PetscErrorCode MatCreateSubMatricesMPI(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
670: 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);}
671: PETSC_EXTERN PetscErrorCode MatDestroyMatrices(PetscInt,Mat *[]);
672: PETSC_EXTERN PetscErrorCode MatDestroySubMatrices(PetscInt,Mat *[]);
673: PETSC_EXTERN PetscErrorCode MatCreateSubMatrix(Mat,IS,IS,MatReuse,Mat *);
674: 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);}
675: PETSC_EXTERN PetscErrorCode MatGetLocalSubMatrix(Mat,IS,IS,Mat*);
676: PETSC_EXTERN PetscErrorCode MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*);
677: PETSC_EXTERN PetscErrorCode MatGetSeqNonzeroStructure(Mat,Mat*);
678: PETSC_EXTERN PetscErrorCode MatDestroySeqNonzeroStructure(Mat*);

680: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJ(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*);
681: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*);
682: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJNumeric(Mat,Mat);
683: PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMat(Mat,MatReuse,Mat*);
684: PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*);
685: PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMatMerge(Mat,MatReuse,IS*,Mat*);
686: PETSC_EXTERN PetscErrorCode MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,Mat*);
687: PETSC_EXTERN PetscErrorCode MatGetGhosts(Mat, PetscInt *,const PetscInt *[]);

689: PETSC_EXTERN PetscErrorCode MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt);
690: PETSC_EXTERN PetscErrorCode MatIncreaseOverlapSplit(Mat mat,PetscInt n,IS is[],PetscInt ov);
691: PETSC_EXTERN PetscErrorCode MatMPIAIJSetUseScalableIncreaseOverlap(Mat,PetscBool);

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

695: PETSC_EXTERN PetscErrorCode MatMatMatMult(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
696: PETSC_EXTERN PetscErrorCode MatGalerkin(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);

698: PETSC_EXTERN PetscErrorCode MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*);
699: PETSC_EXTERN PetscErrorCode MatRARt(Mat,Mat,MatReuse,PetscReal,Mat*);

701: PETSC_EXTERN PetscErrorCode MatTransposeMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
702: PETSC_EXTERN PetscErrorCode MatMatTransposeMult(Mat,Mat,MatReuse,PetscReal,Mat*);

704: PETSC_EXTERN PetscErrorCode MatAXPY(Mat,PetscScalar,Mat,MatStructure);
705: PETSC_EXTERN PetscErrorCode MatAYPX(Mat,PetscScalar,Mat,MatStructure);

707: PETSC_EXTERN PetscErrorCode MatScale(Mat,PetscScalar);
708: PETSC_EXTERN PetscErrorCode MatShift(Mat,PetscScalar);

710: PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
711: PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*);
712: PETSC_EXTERN PetscErrorCode MatGetLayouts(Mat,PetscLayout*,PetscLayout*);
713: PETSC_EXTERN PetscErrorCode MatSetLayouts(Mat,PetscLayout,PetscLayout);
714: PETSC_EXTERN PetscErrorCode MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
715: PETSC_EXTERN PetscErrorCode MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
716: PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
717: PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
718: PETSC_EXTERN PetscErrorCode MatGetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
719: PETSC_EXTERN PetscErrorCode MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
720: PETSC_EXTERN PetscErrorCode MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);

722: PETSC_EXTERN PetscErrorCode MatStashSetInitialSize(Mat,PetscInt,PetscInt);
723: PETSC_EXTERN PetscErrorCode MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);

725: PETSC_EXTERN PetscErrorCode MatInterpolate(Mat,Vec,Vec);
726: PETSC_EXTERN PetscErrorCode MatInterpolateAdd(Mat,Vec,Vec,Vec);
727: PETSC_EXTERN PetscErrorCode MatRestrict(Mat,Vec,Vec);
728: PETSC_EXTERN PetscErrorCode MatMatInterpolate(Mat,Mat,Mat*);
729: PETSC_EXTERN PetscErrorCode MatMatInterpolateAdd(Mat,Mat,Mat,Mat*);
730: PETSC_EXTERN PetscErrorCode MatMatRestrict(Mat,Mat,Mat*);
731: PETSC_EXTERN PetscErrorCode MatCreateVecs(Mat,Vec*,Vec*);
732: 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);}
733: PETSC_EXTERN PetscErrorCode MatCreateRedundantMatrix(Mat,PetscInt,MPI_Comm,MatReuse,Mat*);
734: PETSC_EXTERN PetscErrorCode MatGetMultiProcBlock(Mat,MPI_Comm,MatReuse,Mat*);
735: PETSC_EXTERN PetscErrorCode MatFindZeroDiagonals(Mat,IS*);
736: PETSC_EXTERN PetscErrorCode MatFindOffBlockDiagonalEntries(Mat,IS*);
737: PETSC_EXTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);

739: /*MC
740:    MatSetValue - Set a single entry into a matrix.

742:    Not collective

744:    Synopsis:
745: #include <petscmat.h>
746:      PetscErrorCode MatSetValue(Mat m,PetscInt row,PetscInt col,PetscScalar value,InsertMode mode)

748:    Input Parameters:
749: +  m - the matrix
750: .  row - the row location of the entry
751: .  col - the column location of the entry
752: .  value - the value to insert
753: -  mode - either INSERT_VALUES or ADD_VALUES

755:    Notes:
756:    For efficiency one should use MatSetValues() and set several values simultaneously.

758:    Level: beginner

760: .seealso: MatGetValue(), MatSetValues(), MatSetValueLocal(), MatSetValuesLocal()
761: M*/
762: PETSC_STATIC_INLINE PetscErrorCode MatSetValue(Mat v,PetscInt i,PetscInt j,PetscScalar va,InsertMode mode) {return MatSetValues(v,1,&i,1,&j,&va,mode);}

764: /*@C
765:    MatGetValue - Gets a single value from a matrix

767:    Not Collective; can only return a value owned by the given process

769:    Input Parameters:
770: +  mat - the matrix
771: .  row - the row location of the entry
772: -  col - the column location of the entry

774:    Output Parameter:
775: .  va - the value

777:    Notes:
778:    For efficiency one should use MatGetValues() and get several values simultaneously.

780:    See notes for MatGetValues().

782:    Level: advanced

784: .seealso: MatSetValue(), MatGetValueLocal(), MatGetValues()
785: @*/
786: PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat mat,PetscInt row,PetscInt col,PetscScalar *va) {return MatGetValues(mat,1,&row,1,&col,va);}

788: /*MC
789:    MatSetValueLocal - Inserts or adds a single value into a matrix,
790:    using a local numbering of the nodes.

792:    Not Collective

794:    Input Parameters:
795: +  m - the matrix
796: .  row - the row location of the entry
797: .  col - the column location of the entry
798: .  value - the value to insert
799: -  mode - either INSERT_VALUES or ADD_VALUES

801:    Notes:
802:    For efficiency one should use MatSetValuesLocal() and set several values simultaneously.

804:    See notes for MatSetValuesLocal() for additional information on when and how this function can be used.

806:    Level: intermediate

808: .seealso: MatSetValue(), MatSetValuesLocal()
809: M*/
810: PETSC_STATIC_INLINE PetscErrorCode MatSetValueLocal(Mat v,PetscInt i,PetscInt j,PetscScalar va,InsertMode mode) {return MatSetValuesLocal(v,1,&i,1,&j,&va,mode);}

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

816:    Synopsis:
817: #include <petscmat.h>
818:    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)

820:    Collective

822:    Input Parameters:
823: +  comm - the communicator that will share the eventually allocated matrix
824: .  nrows - the number of LOCAL rows in the matrix
825: -  ncols - the number of LOCAL columns in the matrix

827:    Output Parameters:
828: +  dnz - the array that will be passed to the matrix preallocation routines
829: -  onz - the other array passed to the matrix preallocation routines

831:    Level: intermediate

833:    Notes:
834:     See Users-Manual: ch_performance for more details.

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

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

840: .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
841:           MatPreallocateSymmetricSetLocalBlock()
842: M*/
843: #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
844: do { \
845:   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ncols = (ncols),__rstart,__start,__end; \
846:   _4_PetscCalloc2((size_t)__nrows,&dnz,(size_t)__nrows,&onz);CHKERRQ(_4_ierr); \
847:   __start = 0; __end = __start;                                         \
848:   _4_MPI_Scan(&__ncols,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRMPI(_4_ierr); __start = __end - __ncols;\
849:   _4_MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRMPI(_4_ierr); __rstart = __rstart - __nrows; \
850:   do { } while (0)

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

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

860:    Not Collective

862:    Input Parameters:
863: +  map - the row mapping from local numbering to global numbering
864: .  nrows - the number of rows indicated
865: .  rows - the indices of the rows
866: .  cmap - the column mapping from local to global numbering
867: .  ncols - the number of columns in the matrix
868: .  cols - the columns indicated
869: .  dnz - the array that will be passed to the matrix preallocation routines
870: -  onz - the other array passed to the matrix preallocation routines

872:    Level: intermediate

874:    Notes:
875:     See Users-Manual: ch_performance for more details.

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

879: .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock()
880:           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocalRemoveDups()
881: M*/
882: #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
883: do {\
884:   PetscInt __l;\
885:   _4_ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
886:   _4_ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
887:   for (__l=0;__l<nrows;__l++) {\
888:     _4_MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
889:   }\
890:  } while (0)

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

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

900:    Not Collective

902:    Input Parameters:
903: +  map - the row mapping from local numbering to global numbering
904: .  nrows - the number of rows indicated
905: .  rows - the indices of the rows (these values are mapped to the global values)
906: .  cmap - the column mapping from local to global numbering
907: .  ncols - the number of columns in the matrix   (this value will be changed if duplicate columns are found)
908: .  cols - the columns indicated (these values are mapped to the global values, they are then sorted and duplicates removed)
909: .  dnz - the array that will be passed to the matrix preallocation routines
910: -  onz - the other array passed to the matrix preallocation routines

912:    Level: intermediate

914:    Notes:
915:     See Users-Manual: ch_performance for more details.

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

919: .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock()
920:           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal()
921: M*/
922: #define MatPreallocateSetLocalRemoveDups(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
923: do {\
924:   PetscInt __l;\
925:   _4_ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
926:   _4_ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
927:   _4_PetscSortRemoveDupsInt(&ncols,cols);CHKERRQ(_4_ierr);\
928:   for (__l=0;__l<nrows;__l++) {\
929:     _4_MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
930:   }\
931:  } while (0)

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

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

941:    Not Collective

943:    Input Parameters:
944: +  map - the row mapping from local numbering to global numbering
945: .  nrows - the number of rows indicated
946: .  rows - the indices of the rows
947: .  cmap - the column mapping from local to global numbering
948: .  ncols - the number of columns in the matrix
949: .  cols - the columns indicated
950: .  dnz - the array that will be passed to the matrix preallocation routines
951: -  onz - the other array passed to the matrix preallocation routines

953:    Level: intermediate

955:    Notes:
956:     See Users-Manual: ch_performance for more details.

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

960: .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock()
961:           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock()
962: M*/
963: #define MatPreallocateSetLocalBlock(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
964: do {\
965:   PetscInt __l;\
966:   _4_ISLocalToGlobalMappingApplyBlock(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
967:   _4_ISLocalToGlobalMappingApplyBlock(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
968:   for (__l=0;__l<nrows;__l++) {\
969:     _4_MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
970:   }\
971: } while (0)

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

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

981:    Not Collective

983:    Input Parameters:
984: +  map - the mapping between local numbering and global numbering
985: .  nrows - the number of rows indicated
986: .  rows - the indices of the rows
987: .  ncols - the number of columns in the matrix
988: .  cols - the columns indicated
989: .  dnz - the array that will be passed to the matrix preallocation routines
990: -  onz - the other array passed to the matrix preallocation routines

992:    Level: intermediate

994:    Notes:
995:     See Users-Manual: ch_performance for more details.

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

999: .seealso: MatPreallocateFinalize(), MatPreallocateSet()
1000:           MatPreallocateInitialize(),  MatPreallocateSetLocal()
1001: M*/
1002: #define MatPreallocateSymmetricSetLocalBlock(map,nrows,rows,ncols,cols,dnz,onz) 0;\
1003: do {\
1004:   PetscInt __l;\
1005:   _4_ISLocalToGlobalMappingApplyBlock(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
1006:   _4_ISLocalToGlobalMappingApplyBlock(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
1007:   for (__l=0;__l<nrows;__l++) {\
1008:     _4_MatPreallocateSymmetricSetBlock((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
1009:   }\
1010: } while (0)

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

1016:    Synopsis:
1017: #include <petscmat.h>
1018:    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)

1020:    Not Collective

1022:    Input Parameters:
1023: +  row - the row
1024: .  ncols - the number of columns in the matrix
1025: -  cols - the columns indicated

1027:    Output Parameters:
1028: +  dnz - the array that will be passed to the matrix preallocation routines
1029: -  onz - the other array passed to the matrix preallocation routines

1031:    Level: intermediate

1033:    Notes:
1034:     See Users-Manual: ch_performance for more details.

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

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

1040: .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock()
1041:           MatPreallocateInitialize(), MatPreallocateSetLocal()
1042: M*/
1043: #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
1044: do { PetscInt __i; \
1045:   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);\
1046:   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);\
1047:   for (__i=0; __i<nc; __i++) {\
1048:     if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \
1049:     else if (dnz[row - __rstart] < __ncols) dnz[row - __rstart]++;\
1050:   }\
1051: } while (0)

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

1057:    Synopsis:
1058: #include <petscmat.h>
1059:    PetscErrorCode MatPreallocateSymmetricSetBlock(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)

1061:    Not Collective

1063:    Input Parameters:
1064: +  nrows - the number of rows indicated
1065: .  rows - the indices of the rows
1066: .  ncols - the number of columns in the matrix
1067: .  cols - the columns indicated
1068: .  dnz - the array that will be passed to the matrix preallocation routines
1069: -  onz - the other array passed to the matrix preallocation routines

1071:    Level: intermediate

1073:    Notes:
1074:     See Users-Manual: ch_performance for more details.

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

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

1080: .seealso: MatPreallocateFinalize(), MatPreallocateSet(),  MatPreallocateInitialize(),
1081:           MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal()
1082: M*/
1083: #define MatPreallocateSymmetricSetBlock(row,nc,cols,dnz,onz) 0;\
1084: do { PetscInt __i; \
1085:   for (__i=0; __i<nc; __i++) {\
1086:     if (cols[__i] >= __end) onz[row - __rstart]++; \
1087:     else if (cols[__i] >= row && dnz[row - __rstart] < __ncols) dnz[row - __rstart]++;\
1088:   }\
1089: } while (0)

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

1094:    Synopsis:
1095: #include <petscmat.h>
1096:    PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz)

1098:    Not Collective

1100:    Input Parameters:
1101: +  A - matrix
1102: .  row - row where values exist (must be local to this process)
1103: .  ncols - number of columns
1104: .  cols - columns with nonzeros
1105: .  dnz - the array that will be passed to the matrix preallocation routines
1106: -  onz - the other array passed to the matrix preallocation routines

1108:    Level: intermediate

1110:    Notes:
1111:     See Users-Manual: ch_performance for more details.

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

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

1117: .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
1118:           MatPreallocateSymmetricSetLocalBlock()
1119: M*/
1120: #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)

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

1126:    Synopsis:
1127: #include <petscmat.h>
1128:    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)

1130:    Collective

1132:    Input Parameters:
1133: +  dnz - the array that was be passed to the matrix preallocation routines
1134: -  onz - the other array passed to the matrix preallocation routines

1136:    Level: intermediate

1138:    Notes:
1139:     See Users-Manual: ch_performance for more details.

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

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

1145: .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
1146:           MatPreallocateSymmetricSetLocalBlock()
1147: M*/
1148: #define MatPreallocateFinalize(dnz,onz) 0;_4_PetscFree2(dnz,onz);CHKERRQ(_4_ierr);} while (0)

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

1153: PETSC_EXTERN PetscErrorCode MatInodeAdjustForInodes(Mat,IS*,IS*);
1154: PETSC_EXTERN PetscErrorCode MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);

1156: PETSC_EXTERN PetscErrorCode MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
1157: PETSC_EXTERN PetscErrorCode MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
1158: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
1159: PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
1160: PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
1161: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJFromTriple(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*,PetscInt,PetscBool);

1163: #define MAT_SKIP_ALLOCATION -4

1165: PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
1166: PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
1167: PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
1168: PETSC_EXTERN PetscErrorCode MatSeqAIJSetTotalPreallocation(Mat,PetscInt);

1170: PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1171: PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1172: PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1173: PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []);
1174: PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
1175: PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
1176: PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
1177: PETSC_EXTERN PetscErrorCode MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
1178: PETSC_EXTERN PetscErrorCode MatMPIAdjToSeq(Mat,Mat*);
1179: PETSC_EXTERN PetscErrorCode MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
1180: PETSC_EXTERN PetscErrorCode MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
1181: PETSC_EXTERN PetscErrorCode MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,const PetscInt*[]);
1182: PETSC_EXTERN PetscErrorCode MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,const PetscInt*[]);
1183: PETSC_EXTERN PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat,Mat*);

1185: PETSC_EXTERN PetscErrorCode MatDenseGetLDA(Mat,PetscInt*);
1186: PETSC_EXTERN PetscErrorCode MatDenseSetLDA(Mat,PetscInt);
1187: PETSC_DEPRECATED_FUNCTION("Use MatDenseSetLDA() (since version 3.14)") PETSC_STATIC_INLINE PetscErrorCode MatSeqDenseSetLDA(Mat A,PetscInt lda) {return MatDenseSetLDA(A,lda);}
1188: PETSC_EXTERN PetscErrorCode MatDenseGetLocalMatrix(Mat,Mat*);

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

1192: PETSC_EXTERN PetscErrorCode MatStoreValues(Mat);
1193: PETSC_EXTERN PetscErrorCode MatRetrieveValues(Mat);

1195: PETSC_EXTERN PetscErrorCode MatFindNonzeroRows(Mat,IS*);
1196: PETSC_EXTERN PetscErrorCode MatFindZeroRows(Mat,IS*);
1197: /*
1198:   These routines are not usually accessed directly, rather solving is
1199:   done through the KSP and PC interfaces.
1200: */

1202: /*J
1203:     MatOrderingType - String with the name of a PETSc matrix ordering

1205:    Level: beginner

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

1210: .seealso: MatGetOrdering()
1211: J*/
1212: typedef const char* MatOrderingType;
1213: #define MATORDERINGNATURAL        "natural"
1214: #define MATORDERINGND             "nd"
1215: #define MATORDERING1WD            "1wd"
1216: #define MATORDERINGRCM            "rcm"
1217: #define MATORDERINGQMD            "qmd"
1218: #define MATORDERINGROWLENGTH      "rowlength"
1219: #define MATORDERINGWBM            "wbm"
1220: #define MATORDERINGSPECTRAL       "spectral"
1221: #define MATORDERINGAMD            "amd"            /* only works if UMFPACK is installed with PETSc */
1222: #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 */
1223: #define MATORDERINGEXTERNAL       "external"       /* uses an ordering type internal to the factorization package */

1225: PETSC_EXTERN PetscErrorCode MatGetOrdering(Mat,MatOrderingType,IS*,IS*);
1226: PETSC_EXTERN PetscErrorCode MatGetOrderingList(PetscFunctionList*);
1227: PETSC_EXTERN PetscErrorCode MatOrderingRegister(const char[],PetscErrorCode(*)(Mat,MatOrderingType,IS*,IS*));
1228: PETSC_EXTERN PetscFunctionList MatOrderingList;

1230: PETSC_EXTERN PetscErrorCode MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
1231: PETSC_EXTERN PetscErrorCode MatCreateLaplacian(Mat,PetscReal,PetscBool,Mat*);

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

1235: /*S
1236:     MatFactorShiftType - Numeric Shift for factorizations

1238:    Level: beginner

1240: .seealso: MatGetFactor()
1241: S*/
1242: typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType;
1243: PETSC_EXTERN const char *const MatFactorShiftTypes[];
1244: PETSC_EXTERN const char *const MatFactorShiftTypesDetail[];

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

1249:     Level: beginner

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

1254: .seealso: MatGetFactor()
1255: S*/
1256: typedef enum {MAT_FACTOR_NOERROR,MAT_FACTOR_STRUCT_ZEROPIVOT,MAT_FACTOR_NUMERIC_ZEROPIVOT,MAT_FACTOR_OUTMEMORY,MAT_FACTOR_OTHER} MatFactorError;

1258: PETSC_EXTERN PetscErrorCode MatFactorGetError(Mat,MatFactorError*);
1259: PETSC_EXTERN PetscErrorCode MatFactorClearError(Mat);
1260: PETSC_EXTERN PetscErrorCode MatFactorGetErrorZeroPivot(Mat,PetscReal*,PetscInt*);

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

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

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

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

1273:    Level: developer

1275: .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
1276:           MatFactorInfoInitialize()

1278: S*/
1279: typedef struct {
1280:   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
1281:   PetscReal     usedt;
1282:   PetscReal     dt;             /* drop tolerance */
1283:   PetscReal     dtcol;          /* tolerance for pivoting */
1284:   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
1285:   PetscReal     fill;           /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */
1286:   PetscReal     levels;         /* ICC/ILU(levels) */
1287:   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
1288:                                    factorization may be faster if do not pivot */
1289:   PetscReal     zeropivot;      /* pivot is called zero if less than this */
1290:   PetscReal     shifttype;      /* type of shift added to matrix factor to prevent zero pivots */
1291:   PetscReal     shiftamount;     /* how large the shift is */
1292: } MatFactorInfo;

1294: PETSC_EXTERN PetscErrorCode MatFactorInfoInitialize(MatFactorInfo*);
1295: PETSC_EXTERN PetscErrorCode MatCholeskyFactor(Mat,IS,const MatFactorInfo*);
1296: PETSC_EXTERN PetscErrorCode MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
1297: PETSC_EXTERN PetscErrorCode MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*);
1298: PETSC_EXTERN PetscErrorCode MatLUFactor(Mat,IS,IS,const MatFactorInfo*);
1299: PETSC_EXTERN PetscErrorCode MatILUFactor(Mat,IS,IS,const MatFactorInfo*);
1300: PETSC_EXTERN PetscErrorCode MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
1301: PETSC_EXTERN PetscErrorCode MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
1302: PETSC_EXTERN PetscErrorCode MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
1303: PETSC_EXTERN PetscErrorCode MatICCFactor(Mat,IS,const MatFactorInfo*);
1304: PETSC_EXTERN PetscErrorCode MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*);
1305: PETSC_EXTERN PetscErrorCode MatQRFactor(Mat,IS,const MatFactorInfo*);
1306: PETSC_EXTERN PetscErrorCode MatQRFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
1307: PETSC_EXTERN PetscErrorCode MatQRFactorNumeric(Mat,Mat,const MatFactorInfo*);
1308: PETSC_EXTERN PetscErrorCode MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
1309: PETSC_EXTERN PetscErrorCode MatSolve(Mat,Vec,Vec);
1310: PETSC_EXTERN PetscErrorCode MatForwardSolve(Mat,Vec,Vec);
1311: PETSC_EXTERN PetscErrorCode MatBackwardSolve(Mat,Vec,Vec);
1312: PETSC_EXTERN PetscErrorCode MatSolveAdd(Mat,Vec,Vec,Vec);
1313: PETSC_EXTERN PetscErrorCode MatSolveTranspose(Mat,Vec,Vec);
1314: PETSC_EXTERN PetscErrorCode MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
1315: PETSC_EXTERN PetscErrorCode MatSolves(Mat,Vecs,Vecs);
1316: PETSC_EXTERN PetscErrorCode MatSetUnfactored(Mat);

1318: typedef enum {MAT_FACTOR_SCHUR_UNFACTORED, MAT_FACTOR_SCHUR_FACTORED, MAT_FACTOR_SCHUR_INVERTED} MatFactorSchurStatus;
1319: PETSC_EXTERN PetscErrorCode MatFactorSetSchurIS(Mat,IS);
1320: PETSC_EXTERN PetscErrorCode MatFactorGetSchurComplement(Mat,Mat*,MatFactorSchurStatus*);
1321: PETSC_EXTERN PetscErrorCode MatFactorRestoreSchurComplement(Mat,Mat*,MatFactorSchurStatus);
1322: PETSC_EXTERN PetscErrorCode MatFactorInvertSchurComplement(Mat);
1323: PETSC_EXTERN PetscErrorCode MatFactorCreateSchurComplement(Mat,Mat*,MatFactorSchurStatus*);
1324: PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplement(Mat,Vec,Vec);
1325: PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplementTranspose(Mat,Vec,Vec);
1326: PETSC_EXTERN PetscErrorCode MatFactorFactorizeSchurComplement(Mat);

1328: /*E
1329:     MatSORType - What type of (S)SOR to perform

1331:     Level: beginner

1333:    May be bitwise ORd together

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

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

1339: .seealso: MatSOR()
1340: E*/
1341: typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1342:               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1343:               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
1344:               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
1345: PETSC_EXTERN PetscErrorCode MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);

1347: /*
1348:     These routines are for efficiently computing Jacobians via finite differences.
1349: */

1351: /*S
1352:      MatColoring - Object for managing the coloring of matrices.

1354:    Level: beginner

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

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

1363: .seealso:  MatFDColoringCreate(), MatColoringWeightType, ISColoring, MatFDColoring, DMCreateColoring(), MatColoringCreate(), MatOrdering, MatPartitioning, MatColoringType
1364: S*/
1365: typedef struct _p_MatColoring* MatColoring;

1367: /*J
1368:     MatColoringType - String with the name of a PETSc matrix coloring

1370:    Level: beginner

1372: .seealso: MatColoringSetType(), MatColoring
1373: J*/
1374: typedef const  char*           MatColoringType;
1375: #define MATCOLORINGJP      "jp"
1376: #define MATCOLORINGPOWER   "power"
1377: #define MATCOLORINGNATURAL "natural"
1378: #define MATCOLORINGSL      "sl"
1379: #define MATCOLORINGLF      "lf"
1380: #define MATCOLORINGID      "id"
1381: #define MATCOLORINGGREEDY  "greedy"

1383: /*E
1384:    MatColoringWeightType - Type of weight scheme

1386:     Not Collective

1388: +   MAT_COLORING_RANDOM  - Random weights
1389: .   MAT_COLORING_LEXICAL - Lexical weighting based upon global numbering.
1390: -   MAT_COLORING_LF      - Last-first weighting.

1392:     Level: intermediate

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

1396: .seealso: MatColoring, MatColoringCreate()
1397: E*/
1398: typedef enum {MAT_COLORING_WEIGHT_RANDOM,MAT_COLORING_WEIGHT_LEXICAL,MAT_COLORING_WEIGHT_LF,MAT_COLORING_WEIGHT_SL} MatColoringWeightType;

1400: PETSC_EXTERN PetscErrorCode MatColoringCreate(Mat,MatColoring*);
1401: PETSC_EXTERN PetscErrorCode MatColoringGetDegrees(Mat,PetscInt,PetscInt*);
1402: PETSC_EXTERN PetscErrorCode MatColoringDestroy(MatColoring*);
1403: PETSC_EXTERN PetscErrorCode MatColoringView(MatColoring,PetscViewer);
1404: PETSC_EXTERN PetscErrorCode MatColoringSetType(MatColoring,MatColoringType);
1405: PETSC_EXTERN PetscErrorCode MatColoringSetFromOptions(MatColoring);
1406: PETSC_EXTERN PetscErrorCode MatColoringSetDistance(MatColoring,PetscInt);
1407: PETSC_EXTERN PetscErrorCode MatColoringGetDistance(MatColoring,PetscInt*);
1408: PETSC_EXTERN PetscErrorCode MatColoringSetMaxColors(MatColoring,PetscInt);
1409: PETSC_EXTERN PetscErrorCode MatColoringGetMaxColors(MatColoring,PetscInt*);
1410: PETSC_EXTERN PetscErrorCode MatColoringApply(MatColoring,ISColoring*);
1411: PETSC_EXTERN PetscErrorCode MatColoringRegister(const char[],PetscErrorCode(*)(MatColoring));
1412: PETSC_EXTERN PetscErrorCode MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
1413: PETSC_EXTERN PetscErrorCode MatColoringSetWeightType(MatColoring,MatColoringWeightType);
1414: PETSC_EXTERN PetscErrorCode MatColoringSetWeights(MatColoring,PetscReal*,PetscInt*);
1415: PETSC_EXTERN PetscErrorCode MatColoringCreateWeights(MatColoring,PetscReal **,PetscInt **lperm);
1416: PETSC_EXTERN PetscErrorCode MatColoringTest(MatColoring,ISColoring);
1417: PETSC_DEPRECATED_FUNCTION("Use MatColoringTest() (since version 3.10)") PETSC_STATIC_INLINE PetscErrorCode MatColoringTestValid(MatColoring matcoloring,ISColoring iscoloring) {return MatColoringTest(matcoloring,iscoloring);}
1418: PETSC_EXTERN PetscErrorCode MatISColoringTest(Mat,ISColoring);

1420: /*S
1421:      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1422:         and coloring

1424:    Level: beginner

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

1429: .seealso:  MatFDColoringCreate(), MatColoring, DMCreateColoring()
1430: S*/
1431: typedef struct _p_MatFDColoring* MatFDColoring;

1433: PETSC_EXTERN PetscErrorCode MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
1434: PETSC_EXTERN PetscErrorCode MatFDColoringDestroy(MatFDColoring*);
1435: PETSC_EXTERN PetscErrorCode MatFDColoringView(MatFDColoring,PetscViewer);
1436: PETSC_EXTERN PetscErrorCode MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
1437: PETSC_EXTERN PetscErrorCode MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**);
1438: PETSC_EXTERN PetscErrorCode MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
1439: PETSC_EXTERN PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring);
1440: PETSC_EXTERN PetscErrorCode MatFDColoringApply(Mat,MatFDColoring,Vec,void *);
1441: PETSC_EXTERN PetscErrorCode MatFDColoringSetF(MatFDColoring,Vec);
1442: PETSC_EXTERN PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,const PetscInt*[]);
1443: PETSC_EXTERN PetscErrorCode MatFDColoringSetUp(Mat,ISColoring,MatFDColoring);
1444: PETSC_EXTERN PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring,PetscInt,PetscInt);
1445: PETSC_EXTERN PetscErrorCode MatFDColoringSetValues(Mat,MatFDColoring,const PetscScalar*);

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

1450:    Level: beginner

1452: .seealso:  MatTransposeColoringCreate()
1453: S*/
1454: typedef struct _p_MatTransposeColoring* MatTransposeColoring;

1456: PETSC_EXTERN PetscErrorCode MatTransposeColoringCreate(Mat,ISColoring,MatTransposeColoring *);
1457: PETSC_EXTERN PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring,Mat,Mat);
1458: PETSC_EXTERN PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring,Mat,Mat);
1459: PETSC_EXTERN PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring*);

1461: /*
1462:     These routines are for partitioning matrices: currently used only
1463:   for adjacency matrix, MatCreateMPIAdj().
1464: */

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

1469:    Level: beginner

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

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

1478: .seealso:  MatPartitioningCreate(), MatPartitioningType, MatColoring, MatGetOrdering()
1479: S*/
1480: typedef struct _p_MatPartitioning* MatPartitioning;

1482: /*J
1483:     MatPartitioningType - String with the name of a PETSc matrix partitioning

1485:    Level: beginner
1486: dm
1487: .seealso: MatPartitioningCreate(), MatPartitioning
1488: J*/
1489: typedef const char* MatPartitioningType;
1490: #define MATPARTITIONINGCURRENT  "current"
1491: #define MATPARTITIONINGAVERAGE   "average"
1492: #define MATPARTITIONINGSQUARE   "square"
1493: #define MATPARTITIONINGPARMETIS "parmetis"
1494: #define MATPARTITIONINGCHACO    "chaco"
1495: #define MATPARTITIONINGPARTY    "party"
1496: #define MATPARTITIONINGPTSCOTCH "ptscotch"
1497: #define MATPARTITIONINGHIERARCH  "hierarch"

1499: PETSC_EXTERN PetscErrorCode MatPartitioningCreate(MPI_Comm,MatPartitioning*);
1500: PETSC_EXTERN PetscErrorCode MatPartitioningSetType(MatPartitioning,MatPartitioningType);
1501: PETSC_EXTERN PetscErrorCode MatPartitioningSetNParts(MatPartitioning,PetscInt);
1502: PETSC_EXTERN PetscErrorCode MatPartitioningSetAdjacency(MatPartitioning,Mat);
1503: PETSC_EXTERN PetscErrorCode MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
1504: PETSC_EXTERN PetscErrorCode MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
1505: PETSC_EXTERN PetscErrorCode MatPartitioningSetUseEdgeWeights(MatPartitioning,PetscBool);
1506: PETSC_EXTERN PetscErrorCode MatPartitioningGetUseEdgeWeights(MatPartitioning,PetscBool*);
1507: PETSC_EXTERN PetscErrorCode MatPartitioningApply(MatPartitioning,IS*);
1508: PETSC_EXTERN PetscErrorCode MatPartitioningImprove(MatPartitioning,IS*);
1509: PETSC_EXTERN PetscErrorCode MatPartitioningViewImbalance(MatPartitioning,IS);
1510: PETSC_EXTERN PetscErrorCode MatPartitioningApplyND(MatPartitioning,IS*);
1511: PETSC_EXTERN PetscErrorCode MatPartitioningDestroy(MatPartitioning*);
1512: PETSC_EXTERN PetscErrorCode MatPartitioningRegister(const char[],PetscErrorCode (*)(MatPartitioning));
1513: PETSC_EXTERN PetscErrorCode MatPartitioningView(MatPartitioning,PetscViewer);
1514: PETSC_EXTERN PetscErrorCode MatPartitioningViewFromOptions(MatPartitioning,PetscObject,const char[]);
1515: PETSC_EXTERN PetscErrorCode MatPartitioningSetFromOptions(MatPartitioning);
1516: PETSC_EXTERN PetscErrorCode MatPartitioningGetType(MatPartitioning,MatPartitioningType*);

1518: PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetRepartition(MatPartitioning part);
1519: PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
1520: PETSC_EXTERN PetscErrorCode MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *);

1522: typedef enum { MP_CHACO_MULTILEVEL=1,MP_CHACO_SPECTRAL=2,MP_CHACO_LINEAR=4,MP_CHACO_RANDOM=5,MP_CHACO_SCATTERED=6 } MPChacoGlobalType;
1523: PETSC_EXTERN const char *const MPChacoGlobalTypes[];
1524: typedef enum { MP_CHACO_KERNIGHAN=1,MP_CHACO_NONE=2 } MPChacoLocalType;
1525: PETSC_EXTERN const char *const MPChacoLocalTypes[];
1526: typedef enum { MP_CHACO_LANCZOS=0,MP_CHACO_RQI=1 } MPChacoEigenType;
1527: PETSC_EXTERN const char *const MPChacoEigenTypes[];

1529: PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetGlobal(MatPartitioning,MPChacoGlobalType);
1530: PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetGlobal(MatPartitioning,MPChacoGlobalType*);
1531: PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetLocal(MatPartitioning,MPChacoLocalType);
1532: PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetLocal(MatPartitioning,MPChacoLocalType*);
1533: PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
1534: PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
1535: PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenSolver(MatPartitioning,MPChacoEigenType*);
1536: PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenTol(MatPartitioning,PetscReal);
1537: PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenTol(MatPartitioning,PetscReal*);
1538: PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenNumber(MatPartitioning,PetscInt);
1539: PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenNumber(MatPartitioning,PetscInt*);

1541: #define MP_PARTY_OPT "opt"
1542: #define MP_PARTY_LIN "lin"
1543: #define MP_PARTY_SCA "sca"
1544: #define MP_PARTY_RAN "ran"
1545: #define MP_PARTY_GBF "gbf"
1546: #define MP_PARTY_GCF "gcf"
1547: #define MP_PARTY_BUB "bub"
1548: #define MP_PARTY_DEF "def"
1549: PETSC_EXTERN PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning,const char*);
1550: #define MP_PARTY_HELPFUL_SETS "hs"
1551: #define MP_PARTY_KERNIGHAN_LIN "kl"
1552: #define MP_PARTY_NONE "no"
1553: PETSC_EXTERN PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning,const char*);
1554: PETSC_EXTERN PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
1555: PETSC_EXTERN PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning,PetscBool);
1556: PETSC_EXTERN PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool);

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

1561: PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning,PetscReal);
1562: PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning,PetscReal*);
1563: PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning,MPPTScotchStrategyType);
1564: PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning,MPPTScotchStrategyType*);

1566: /*
1567:  * hierarchical partitioning
1568:  */
1569: PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetFineparts(MatPartitioning,IS*);
1570: PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetCoarseparts(MatPartitioning,IS*);
1571: PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNcoarseparts(MatPartitioning,PetscInt);
1572: PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNfineparts(MatPartitioning, PetscInt);

1574: PETSC_EXTERN PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*);
1575: PETSC_EXTERN PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*);

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

1744: /*
1745:    Codes for matrices stored on disk. By default they are
1746:    stored in a universal format. By changing the format with
1747:    PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will
1748:    be stored in a way natural for the matrix, for example dense matrices
1749:    would be stored as dense. Matrices stored this way may only be
1750:    read into matrices of the same type.
1751: */
1752: #define MATRIX_BINARY_FORMAT_DENSE -1

1754: PETSC_EXTERN PetscErrorCode MatMPIBAIJSetHashTableFactor(Mat,PetscReal);

1756: PETSC_EXTERN PetscErrorCode MatISSetLocalMatType(Mat,MatType);
1757: PETSC_EXTERN PetscErrorCode MatISSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1758: PETSC_EXTERN PetscErrorCode MatISStoreL2L(Mat,PetscBool);
1759: PETSC_EXTERN PetscErrorCode MatISFixLocalEmpty(Mat,PetscBool);
1760: PETSC_EXTERN PetscErrorCode MatISGetLocalMat(Mat,Mat*);
1761: PETSC_EXTERN PetscErrorCode MatISRestoreLocalMat(Mat,Mat*);
1762: PETSC_EXTERN PetscErrorCode MatISSetLocalMat(Mat,Mat);
1763: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use the MatConvert() interface (since version 3.10)") PetscErrorCode MatISGetMPIXAIJ(Mat,MatReuse,Mat*);

1765: /*S
1766:      MatNullSpace - Object that removes a null space from a vector, i.e.
1767:          orthogonalizes the vector to a subspace

1769:    Level: advanced

1771:   Users manual sections:
1772: .   sec_singular

1774: .seealso:  MatNullSpaceCreate()
1775: S*/
1776: typedef struct _p_MatNullSpace* MatNullSpace;

1778: PETSC_EXTERN PetscErrorCode MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*);
1779: PETSC_EXTERN PetscErrorCode MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*);
1780: PETSC_EXTERN PetscErrorCode MatNullSpaceDestroy(MatNullSpace*);
1781: PETSC_EXTERN PetscErrorCode MatNullSpaceRemove(MatNullSpace,Vec);
1782: PETSC_EXTERN PetscErrorCode MatGetNullSpace(Mat, MatNullSpace *);
1783: PETSC_EXTERN PetscErrorCode MatGetTransposeNullSpace(Mat, MatNullSpace *);
1784: PETSC_EXTERN PetscErrorCode MatSetTransposeNullSpace(Mat,MatNullSpace);
1785: PETSC_EXTERN PetscErrorCode MatSetNullSpace(Mat,MatNullSpace);
1786: PETSC_EXTERN PetscErrorCode MatSetNearNullSpace(Mat,MatNullSpace);
1787: PETSC_EXTERN PetscErrorCode MatGetNearNullSpace(Mat,MatNullSpace*);
1788: PETSC_EXTERN PetscErrorCode MatNullSpaceTest(MatNullSpace,Mat,PetscBool  *);
1789: PETSC_EXTERN PetscErrorCode MatNullSpaceView(MatNullSpace,PetscViewer);
1790: PETSC_EXTERN PetscErrorCode MatNullSpaceGetVecs(MatNullSpace,PetscBool*,PetscInt*,const Vec**);
1791: PETSC_EXTERN PetscErrorCode MatNullSpaceCreateRigidBody(Vec,MatNullSpace*);

1793: PETSC_EXTERN PetscErrorCode MatReorderingSeqSBAIJ(Mat,IS);
1794: PETSC_EXTERN PetscErrorCode MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
1795: PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);

1797: PETSC_EXTERN PetscErrorCode MatCreateMAIJ(Mat,PetscInt,Mat*);
1798: PETSC_EXTERN PetscErrorCode MatMAIJRedimension(Mat,PetscInt,Mat*);
1799: PETSC_EXTERN PetscErrorCode MatMAIJGetAIJ(Mat,Mat*);

1801: PETSC_EXTERN PetscErrorCode MatComputeOperator(Mat,MatType,Mat*);
1802: PETSC_EXTERN PetscErrorCode MatComputeOperatorTranspose(Mat,MatType,Mat*);

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

1807: PETSC_EXTERN PetscErrorCode MatCreateKAIJ(Mat,PetscInt,PetscInt,const PetscScalar[],const PetscScalar[],Mat*);
1808: PETSC_EXTERN PetscErrorCode MatKAIJGetAIJ(Mat,Mat*);
1809: PETSC_EXTERN PetscErrorCode MatKAIJGetS(Mat,PetscInt*,PetscInt*,PetscScalar**);
1810: PETSC_EXTERN PetscErrorCode MatKAIJGetSRead(Mat,PetscInt*,PetscInt*,const PetscScalar**);
1811: PETSC_EXTERN PetscErrorCode MatKAIJRestoreS(Mat,PetscScalar**);
1812: PETSC_EXTERN PetscErrorCode MatKAIJRestoreSRead(Mat,const PetscScalar**);
1813: PETSC_EXTERN PetscErrorCode MatKAIJGetT(Mat,PetscInt*,PetscInt*,PetscScalar**);
1814: PETSC_EXTERN PetscErrorCode MatKAIJGetTRead(Mat,PetscInt*,PetscInt*,const PetscScalar**);
1815: PETSC_EXTERN PetscErrorCode MatKAIJRestoreT(Mat,PetscScalar**);
1816: PETSC_EXTERN PetscErrorCode MatKAIJRestoreTRead(Mat,const PetscScalar**);
1817: PETSC_EXTERN PetscErrorCode MatKAIJSetAIJ(Mat,Mat);
1818: PETSC_EXTERN PetscErrorCode MatKAIJSetS(Mat,PetscInt,PetscInt,const PetscScalar[]);
1819: PETSC_EXTERN PetscErrorCode MatKAIJSetT(Mat,PetscInt,PetscInt,const PetscScalar[]);
1820: PETSC_EXTERN PetscErrorCode MatKAIJGetScaledIdentity(Mat,PetscBool*);

1822: PETSC_EXTERN PetscErrorCode MatDiagonalScaleLocal(Mat,Vec);

1824: PETSC_EXTERN PetscErrorCode MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
1825: PETSC_EXTERN PetscErrorCode MatMFFDSetBase(Mat,Vec,Vec);
1826: PETSC_EXTERN PetscErrorCode MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*);
1827: PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*));
1828: PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec));
1829: PETSC_EXTERN PetscErrorCode MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt);
1830: PETSC_EXTERN PetscErrorCode MatMFFDResetHHistory(Mat);
1831: PETSC_EXTERN PetscErrorCode MatMFFDSetFunctionError(Mat,PetscReal);
1832: PETSC_EXTERN PetscErrorCode MatMFFDSetPeriod(Mat,PetscInt);
1833: PETSC_EXTERN PetscErrorCode MatMFFDGetH(Mat,PetscScalar *);
1834: PETSC_EXTERN PetscErrorCode MatMFFDSetOptionsPrefix(Mat,const char[]);
1835: PETSC_EXTERN PetscErrorCode MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*);
1836: PETSC_EXTERN PetscErrorCode MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*);

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

1842:     Notes:
1843:     MATMFFD is a specific MatType which uses the MatMFFD data structure

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

1847:     Level: developer

1849: .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister()
1850: S*/
1851: typedef struct _p_MatMFFD* MatMFFD;

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

1856:    Level: beginner

1858: .seealso: MatMFFDSetType(), MatMFFDRegister()
1859: J*/
1860: typedef const char* MatMFFDType;
1861: #define MATMFFD_DS  "ds"
1862: #define MATMFFD_WP  "wp"

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

1867: PETSC_EXTERN PetscErrorCode MatMFFDDSSetUmin(Mat,PetscReal);
1868: PETSC_EXTERN PetscErrorCode MatMFFDWPSetComputeNormU(Mat,PetscBool);

1870: PETSC_EXTERN PetscErrorCode MatFDColoringSetType(MatFDColoring,MatMFFDType);

1872: PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
1873: PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);

1875: #ifdef PETSC_HAVE_H2OPUS
1876: PETSC_EXTERN_TYPEDEF typedef PetscScalar (*MatH2OpusKernel)(PetscInt,PetscReal[],PetscReal[],void*);
1877: PETSC_EXTERN PetscErrorCode MatCreateH2OpusFromKernel(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscReal[],PetscBool,MatH2OpusKernel,void*,PetscReal,PetscInt,PetscInt,Mat*);
1878: PETSC_EXTERN PetscErrorCode MatCreateH2OpusFromMat(Mat,PetscInt,const PetscReal[],PetscBool,PetscReal,PetscInt,PetscInt,PetscInt,PetscReal,Mat*);
1879: PETSC_EXTERN PetscErrorCode MatH2OpusSetSamplingMat(Mat,Mat,PetscInt,PetscReal);
1880: PETSC_EXTERN PetscErrorCode MatH2OpusOrthogonalize(Mat);
1881: PETSC_EXTERN PetscErrorCode MatH2OpusCompress(Mat,PetscReal);
1882: PETSC_EXTERN PetscErrorCode MatH2OpusSetNativeMult(Mat,PetscBool);
1883: PETSC_EXTERN PetscErrorCode MatH2OpusGetNativeMult(Mat,PetscBool*);
1884: PETSC_EXTERN PetscErrorCode MatH2OpusGetIndexMap(Mat,IS*);
1885: PETSC_EXTERN PetscErrorCode MatH2OpusMapVec(Mat,PetscBool,Vec,Vec*);
1886: #endif

1888: #ifdef PETSC_HAVE_HTOOL
1889: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*MatHtoolKernel)(PetscInt,PetscInt,PetscInt,const PetscInt*,const PetscInt*,PetscScalar*,void*);
1890: PETSC_EXTERN PetscErrorCode MatCreateHtoolFromKernel(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscReal[],const PetscReal[],MatHtoolKernel,void*,Mat*);
1891: PETSC_EXTERN PetscErrorCode MatHtoolSetKernel(Mat,MatHtoolKernel,void*);
1892: PETSC_EXTERN PetscErrorCode MatHtoolGetPermutationSource(Mat,IS*);
1893: PETSC_EXTERN PetscErrorCode MatHtoolGetPermutationTarget(Mat,IS*);
1894: PETSC_EXTERN PetscErrorCode MatHtoolUsePermutation(Mat,PetscBool);
1895: /*E
1896:      MatHtoolCompressorType - Indicates the type of compressor used by a MATHTOOL

1898:    Level: beginner

1900:     Values:
1901: +   MAT_HTOOL_COMPRESSOR_SYMPARTIAL_ACA (default) - symmetric partial adaptive cross approximation
1902: .   MAT_HTOOL_COMPRESSOR_FULL_ACA - full adaptive cross approximation
1903: -   MAT_HTOOL_COMPRESSOR_SVD - singular value decomposition

1905: .seealso: MatCreateHtoolFromKernel(), MATHTOOL, MatHtoolClusteringType
1906: E*/
1907: typedef enum { MAT_HTOOL_COMPRESSOR_SYMPARTIAL_ACA, MAT_HTOOL_COMPRESSOR_FULL_ACA, MAT_HTOOL_COMPRESSOR_SVD } MatHtoolCompressorType;
1908: /*E
1909:      MatHtoolClusteringType - Indicates the type of clustering used by a MATHTOOL

1911:    Level: beginner

1913:     Values:
1914: +   MAT_HTOOL_CLUSTERING_PCA_REGULAR (default) - axis computed via principle component analysis, split uniformly
1915: .   MAT_HTOOL_CLUSTERING_PCA_GEOMETRIC - axis computed via principle component analysis, split barycentrically
1916: .   MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_REGULAR - axis along the largest extent of the bounding box, split uniformly
1917: -   MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_GEOMETRIC - axis along the largest extent of the bounding box, split barycentrically

1919:     Notes: higher-dimensional clustering is not yet supported in Htool, but once it is, one should add BOUNDING_BOX_{2,3} types

1921: .seealso: MatCreateHtoolFromKernel(), MATHTOOL, MatHtoolCompressorType
1922: E*/
1923: typedef enum { MAT_HTOOL_CLUSTERING_PCA_REGULAR, MAT_HTOOL_CLUSTERING_PCA_GEOMETRIC, MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_REGULAR, MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_GEOMETRIC } MatHtoolClusteringType;
1924: #endif

1926: /*
1927:    PETSc interface to MUMPS
1928: */
1929: #ifdef PETSC_HAVE_MUMPS
1930: PETSC_EXTERN PetscErrorCode MatMumpsSetIcntl(Mat,PetscInt,PetscInt);
1931: PETSC_EXTERN PetscErrorCode MatMumpsGetIcntl(Mat,PetscInt,PetscInt*);
1932: PETSC_EXTERN PetscErrorCode MatMumpsSetCntl(Mat,PetscInt,PetscReal);
1933: PETSC_EXTERN PetscErrorCode MatMumpsGetCntl(Mat,PetscInt,PetscReal*);

1935: PETSC_EXTERN PetscErrorCode MatMumpsGetInfo(Mat,PetscInt,PetscInt*);
1936: PETSC_EXTERN PetscErrorCode MatMumpsGetInfog(Mat,PetscInt,PetscInt*);
1937: PETSC_EXTERN PetscErrorCode MatMumpsGetRinfo(Mat,PetscInt,PetscReal*);
1938: PETSC_EXTERN PetscErrorCode MatMumpsGetRinfog(Mat,PetscInt,PetscReal*);
1939: PETSC_EXTERN PetscErrorCode MatMumpsGetInverse(Mat,Mat);
1940: PETSC_EXTERN PetscErrorCode MatMumpsGetInverseTranspose(Mat,Mat);
1941: #endif

1943: /*
1944:    PETSc interface to Mkl_Pardiso
1945: */
1946: #ifdef PETSC_HAVE_MKL_PARDISO
1947: PETSC_EXTERN PetscErrorCode MatMkl_PardisoSetCntl(Mat,PetscInt,PetscInt);
1948: #endif

1950: /*
1951:    PETSc interface to Mkl_CPardiso
1952: */
1953: #ifdef PETSC_HAVE_MKL_CPARDISO
1954: PETSC_EXTERN PetscErrorCode MatMkl_CPardisoSetCntl(Mat,PetscInt,PetscInt);
1955: #endif

1957: /*
1958:    PETSc interface to SUPERLU
1959: */
1960: #ifdef PETSC_HAVE_SUPERLU
1961: PETSC_EXTERN PetscErrorCode MatSuperluSetILUDropTol(Mat,PetscReal);
1962: #endif

1964: /*
1965:    PETSc interface to SUPERLU_DIST
1966: */
1967: #ifdef PETSC_HAVE_SUPERLU_DIST
1968: PETSC_EXTERN PetscErrorCode MatSuperluDistGetDiagU(Mat,PetscScalar*);
1969: #endif

1971: /*
1972:    PETSc interface to STRUMPACK
1973: */
1974: #ifdef PETSC_HAVE_STRUMPACK
1975: /*E
1976:     MatSTRUMPACKReordering - sparsity reducing ordering to be used in STRUMPACK

1978:     Level: intermediate
1979: E*/
1980: typedef enum {MAT_STRUMPACK_NATURAL=0,
1981:               MAT_STRUMPACK_METIS=1,
1982:               MAT_STRUMPACK_PARMETIS=2,
1983:               MAT_STRUMPACK_SCOTCH=3,
1984:               MAT_STRUMPACK_PTSCOTCH=4,
1985:               MAT_STRUMPACK_RCM=5} MatSTRUMPACKReordering;
1986: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetReordering(Mat,MatSTRUMPACKReordering);
1987: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetColPerm(Mat,PetscBool);
1988: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSRelTol(Mat,PetscReal);
1989: PETSC_DEPRECATED_FUNCTION("Use MatSTRUMPACKSetHSSRelTol() (since version 3.9)") PETSC_STATIC_INLINE PetscErrorCode MatSTRUMPACKSetHSSRelCompTol(Mat mat,PetscReal rtol) {return MatSTRUMPACKSetHSSRelTol(mat,rtol);}
1990: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSAbsTol(Mat,PetscReal);
1991: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSMinSepSize(Mat,PetscInt);
1992: PETSC_DEPRECATED_FUNCTION("Use MatSTRUMPACKSetHSSMinSepSize() (since version 3.9)") PETSC_STATIC_INLINE PetscErrorCode MatSTRUMPACKSetHSSMinSize(Mat mat,PetscInt hssminsize) {return MatSTRUMPACKSetHSSMinSepSize(mat,hssminsize);}
1993: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSMaxRank(Mat,PetscInt);
1994: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSLeafSize(Mat,PetscInt);
1995: #endif

1997: PETSC_EXTERN PetscErrorCode MatBindToCPU(Mat,PetscBool);
1998: PETSC_EXTERN PetscErrorCode MatBoundToCPU(Mat,PetscBool*);
1999: PETSC_DEPRECATED_FUNCTION("Use MatBindToCPU (since v3.13)") PETSC_STATIC_INLINE PetscErrorCode MatPinToCPU(Mat A,PetscBool flg) {return MatBindToCPU(A,flg);}

2001: typedef struct _n_SplitCSRMat *PetscSplitCSRDataStructure;
2002: PETSC_EXTERN PetscErrorCode MatCUSPARSEGetDeviceMatWrite(Mat,PetscSplitCSRDataStructure*);

2004: #ifdef PETSC_HAVE_KOKKOS_KERNELS
2005: PETSC_EXTERN PetscErrorCode MatKokkosGetDeviceMatWrite(Mat,PetscSplitCSRDataStructure*);
2006: PETSC_EXTERN PetscErrorCode MatSeqAIJKokkosSetDeviceMat(Mat,PetscSplitCSRDataStructure);
2007: PETSC_EXTERN PetscErrorCode MatSeqAIJKokkosGetDeviceMat(Mat,PetscSplitCSRDataStructure*);
2008: #endif

2010: #ifdef PETSC_HAVE_CUDA
2011: /*E
2012:     MatCUSPARSEStorageFormat - indicates the storage format for CUSPARSE (GPU)
2013:     matrices.

2015:     Not Collective

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

2021:     Level: intermediate

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

2025: .seealso: MatCUSPARSESetFormat(), MatCUSPARSEFormatOperation
2026: E*/

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

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

2033: /*E
2034:     MatCUSPARSEFormatOperation - indicates the operation of CUSPARSE (GPU)
2035:     matrices whose operation should use a particular storage format.

2037:     Not Collective

2039: +   MAT_CUSPARSE_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult
2040: .   MAT_CUSPARSE_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult
2041: .   MAT_CUSPARSE_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult
2042: -   MAT_CUSPARSE_ALL - sets the storage format for all CUSPARSE (GPU) matrices

2044:     Level: intermediate

2046: .seealso: MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat
2047: E*/
2048: typedef enum {MAT_CUSPARSE_MULT_DIAG, MAT_CUSPARSE_MULT_OFFDIAG, MAT_CUSPARSE_MULT, MAT_CUSPARSE_ALL} MatCUSPARSEFormatOperation;

2050: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
2051: PETSC_EXTERN PetscErrorCode MatCreateAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
2052: PETSC_EXTERN PetscErrorCode MatCUSPARSESetFormat(Mat,MatCUSPARSEFormatOperation,MatCUSPARSEStorageFormat);
2053: typedef struct Mat_SeqAIJCUSPARSETriFactors* Mat_SeqAIJCUSPARSETriFactors_p;
2054: PETSC_EXTERN PetscErrorCode MatSeqAIJCUSPARSEGetIJ(Mat,PetscBool,const int**,const int**);
2055: PETSC_EXTERN PetscErrorCode MatSeqAIJCUSPARSERestoreIJ(Mat,PetscBool,const int**,const int**);
2056: PETSC_EXTERN PetscErrorCode MatSeqAIJCUSPARSEGetArrayRead(Mat,const PetscScalar**);
2057: PETSC_EXTERN PetscErrorCode MatSeqAIJCUSPARSERestoreArrayRead(Mat,const PetscScalar**);
2058: PETSC_EXTERN PetscErrorCode MatSeqAIJCUSPARSEGetArrayWrite(Mat,PetscScalar**);
2059: PETSC_EXTERN PetscErrorCode MatSeqAIJCUSPARSERestoreArrayWrite(Mat,PetscScalar**);
2060: PETSC_EXTERN PetscErrorCode MatSeqAIJCUSPARSEGetArray(Mat,PetscScalar**);
2061: PETSC_EXTERN PetscErrorCode MatSeqAIJCUSPARSERestoreArray(Mat,PetscScalar**);

2063: PETSC_EXTERN PetscErrorCode MatCreateDenseCUDA(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
2064: PETSC_EXTERN PetscErrorCode MatCreateSeqDenseCUDA(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
2065: PETSC_EXTERN PetscErrorCode MatMPIDenseCUDASetPreallocation(Mat,PetscScalar[]);
2066: PETSC_EXTERN PetscErrorCode MatSeqDenseCUDASetPreallocation(Mat,PetscScalar[]);
2067: PETSC_EXTERN PetscErrorCode MatDenseCUDAGetArrayWrite(Mat,PetscScalar**);
2068: PETSC_EXTERN PetscErrorCode MatDenseCUDAGetArrayRead(Mat,const PetscScalar**);
2069: PETSC_EXTERN PetscErrorCode MatDenseCUDAGetArray(Mat,PetscScalar**);
2070: PETSC_EXTERN PetscErrorCode MatDenseCUDARestoreArrayWrite(Mat,PetscScalar**);
2071: PETSC_EXTERN PetscErrorCode MatDenseCUDARestoreArrayRead(Mat,const PetscScalar**);
2072: PETSC_EXTERN PetscErrorCode MatDenseCUDARestoreArray(Mat,PetscScalar**);
2073: PETSC_EXTERN PetscErrorCode MatDenseCUDAPlaceArray(Mat,const PetscScalar*);
2074: PETSC_EXTERN PetscErrorCode MatDenseCUDAReplaceArray(Mat,const PetscScalar*);
2075: PETSC_EXTERN PetscErrorCode MatDenseCUDAResetArray(Mat);

2077: #endif

2079: #if defined(PETSC_HAVE_VIENNACL)
2080: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
2081: PETSC_EXTERN PetscErrorCode MatCreateAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
2082: #endif

2084: /*
2085:    PETSc interface to FFTW
2086: */
2087: #if defined(PETSC_HAVE_FFTW)
2088: PETSC_EXTERN PetscErrorCode VecScatterPetscToFFTW(Mat,Vec,Vec);
2089: PETSC_EXTERN PetscErrorCode VecScatterFFTWToPetsc(Mat,Vec,Vec);
2090: PETSC_EXTERN PetscErrorCode MatCreateVecsFFTW(Mat,Vec*,Vec*,Vec*);
2091: #endif

2093: #if defined(PETSC_HAVE_SCALAPACK)
2094: PETSC_EXTERN PetscErrorCode MatCreateScaLAPACK(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
2095: PETSC_EXTERN PetscErrorCode MatScaLAPACKSetBlockSizes(Mat,PetscInt,PetscInt);
2096: PETSC_EXTERN PetscErrorCode MatScaLAPACKGetBlockSizes(Mat,PetscInt*,PetscInt*);
2097: #endif

2099: PETSC_EXTERN PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*);
2100: PETSC_EXTERN PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*);
2101: PETSC_EXTERN PetscErrorCode MatNestGetISs(Mat,IS[],IS[]);
2102: PETSC_EXTERN PetscErrorCode MatNestGetLocalISs(Mat,IS[],IS[]);
2103: PETSC_EXTERN PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***);
2104: PETSC_EXTERN PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*);
2105: PETSC_EXTERN PetscErrorCode MatNestSetVecType(Mat,VecType);
2106: PETSC_EXTERN PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]);
2107: PETSC_EXTERN PetscErrorCode MatNestSetSubMat(Mat,PetscInt,PetscInt,Mat);

2109: PETSC_EXTERN PetscErrorCode MatChop(Mat,PetscReal);
2110: PETSC_EXTERN PetscErrorCode MatComputeBandwidth(Mat,PetscReal,PetscInt*);

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

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

2116: PETSC_INTERN PetscErrorCode MatHeaderMerge(Mat,Mat*);
2117: PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat,Mat*);

2119: #endif