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)") 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)") 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)") 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: 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, MatProductSetAlgorithm(), MatProductType
196: J*/
197: typedef const char* MatProductAlgorithm;
198: #define MATPRODUCTALGORITHMDEFAULT "default"
199: #define MATPRODUCTALGORITHMSORTED "sorted"
200: #define MATPRODUCTALGORITHMSCALABLE "scalable"
201: #define MATPRODUCTALGORITHMSCALABLEFAST "scalable_fast"
202: #define MATPRODUCTALGORITHMHEAP "heap"
203: #define MATPRODUCTALGORITHMBHEAP "btheap"
204: #define MATPRODUCTALGORITHMLLCONDENSED "llcondensed"
205: #define MATPRODUCTALGORITHMROWMERGE "rowmerge"
206: #define MATPRODUCTALGORITHMOUTERPRODUCT "outerproduct"
207: #define MATPRODUCTALGORITHMATB "at*b"
208: #define MATPRODUCTALGORITHMRAP "rap"
209: #define MATPRODUCTALGORITHMNONSCALABLE "nonscalable"
210: #define MATPRODUCTALGORITHMSEQMPI "seqmpi"
211: #define MATPRODUCTALGORITHMBACKEND "backend"
212: #define MATPRODUCTALGORITHMOVERLAPPING "overlapping"
213: #define MATPRODUCTALGORITHMMERGED "merged"
214: #define MATPRODUCTALGORITHMALLATONCE "allatonce"
215: #define MATPRODUCTALGORITHMALLATONCEMERGED "allatonce_merged"
216: #define MATPRODUCTALGORITHMALLGATHERV "allgatherv"
217: #define MATPRODUCTALGORITHMCYCLIC "cyclic"
218: #if defined(PETSC_HAVE_HYPRE)
219: #define MATPRODUCTALGORITHMHYPRE "hypre"
220: #endif

222: PETSC_EXTERN PetscErrorCode MatProductCreate(Mat,Mat,Mat,Mat*);
223: PETSC_EXTERN PetscErrorCode MatProductCreateWithMat(Mat,Mat,Mat,Mat);
224: PETSC_EXTERN PetscErrorCode MatProductSetType(Mat,MatProductType);
225: PETSC_EXTERN PetscErrorCode MatProductSetAlgorithm(Mat,MatProductAlgorithm);
226: PETSC_EXTERN PetscErrorCode MatProductSetFill(Mat,PetscReal);
227: PETSC_EXTERN PetscErrorCode MatProductSetFromOptions(Mat);
228: PETSC_EXTERN PetscErrorCode MatProductSymbolic(Mat);
229: PETSC_EXTERN PetscErrorCode MatProductNumeric(Mat);
230: PETSC_EXTERN PetscErrorCode MatProductReplaceMats(Mat,Mat,Mat,Mat);
231: PETSC_EXTERN PetscErrorCode MatProductClear(Mat);
232: PETSC_EXTERN PetscErrorCode MatProductView(Mat,PetscViewer);
233: PETSC_EXTERN PetscErrorCode MatProductGetType(Mat,MatProductType*);
234: PETSC_EXTERN PetscErrorCode MatProductGetMats(Mat,Mat*,Mat*,Mat*);

236: /* Logging support */
237: #define    MAT_FILE_CLASSID 1211216    /* used to indicate matrices in binary files */
238: PETSC_EXTERN PetscClassId MAT_CLASSID;
239: PETSC_EXTERN PetscClassId MAT_COLORING_CLASSID;
240: PETSC_EXTERN PetscClassId MAT_FDCOLORING_CLASSID;
241: PETSC_EXTERN PetscClassId MAT_TRANSPOSECOLORING_CLASSID;
242: PETSC_EXTERN PetscClassId MAT_PARTITIONING_CLASSID;
243: PETSC_EXTERN PetscClassId MAT_COARSEN_CLASSID;
244: PETSC_EXTERN PetscClassId MAT_NULLSPACE_CLASSID;
245: PETSC_EXTERN PetscClassId MATMFFD_CLASSID;

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

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

256:     Level: beginner

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

260: .seealso: MatCreateSubMatrices(), MatCreateSubMatrix(), MatDestroyMatrices(), MatConvert()
261: E*/
262: typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX,MAT_IGNORE_MATRIX,MAT_INPLACE_MATRIX} MatReuse;

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

268:     Level: beginner

270: .seealso: MatGetSeqNonzeroStructure()
271: E*/
272: typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatCreateSubMatrixOption;

274: PETSC_EXTERN PetscErrorCode MatInitializePackage(void);

276: PETSC_EXTERN PetscErrorCode MatCreate(MPI_Comm,Mat*);
277: PETSC_EXTERN PetscErrorCode MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt);
278: PETSC_EXTERN PetscErrorCode MatSetType(Mat,MatType);
279: PETSC_EXTERN PetscErrorCode MatGetVecType(Mat,VecType*);
280: PETSC_EXTERN PetscErrorCode MatSetVecType(Mat,VecType);
281: PETSC_EXTERN PetscErrorCode MatSetFromOptions(Mat);
282: PETSC_EXTERN PetscErrorCode MatViewFromOptions(Mat,PetscObject,const char[]);
283: PETSC_EXTERN PetscErrorCode MatRegister(const char[],PetscErrorCode(*)(Mat));
284: PETSC_EXTERN PetscErrorCode MatRegisterRootName(const char[],const char[],const char[]);
285: PETSC_EXTERN PetscErrorCode MatSetOptionsPrefix(Mat,const char[]);
286: PETSC_EXTERN PetscErrorCode MatAppendOptionsPrefix(Mat,const char[]);
287: PETSC_EXTERN PetscErrorCode MatGetOptionsPrefix(Mat,const char*[]);
288: PETSC_EXTERN PetscErrorCode MatSetErrorIfFailure(Mat,PetscBool);

290: PETSC_EXTERN PetscFunctionList MatList;
291: PETSC_EXTERN PetscFunctionList MatColoringList;
292: PETSC_EXTERN PetscFunctionList MatPartitioningList;

294: /*E
295:     MatStructure - Indicates if two matrices have the same nonzero structure

297:     Level: beginner

299: $  SAME_NONZERO_PATTERN  - the two matrices have identical nonzero patterns
300: $  DIFFERENT_NONZERO_PATTERN - the two matrices may have different nonzero patterns
301: $  SUBSET_NONZERO_PATTERN - the nonzero pattern of the second matrix is a subset of the nonzero pattern of the first matrix
302: $  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

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

307: .seealso: MatCopy(), MatAXPY(), MatAYPX()
308: E*/
309: typedef enum {DIFFERENT_NONZERO_PATTERN,SUBSET_NONZERO_PATTERN,SAME_NONZERO_PATTERN,UNKNOWN_NONZERO_PATTERN} MatStructure;
310: PETSC_EXTERN const char *const MatStructures[];

312: #if defined PETSC_HAVE_MKL_SPARSE
313: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJMKL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
314: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJMKL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
315: PETSC_EXTERN PetscErrorCode MatCreateBAIJMKL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
316: PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJMKL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
317: #endif

319: PETSC_EXTERN PetscErrorCode MatCreateSeqSELL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
320: PETSC_EXTERN PetscErrorCode MatCreateSELL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
321: PETSC_EXTERN PetscErrorCode MatSeqSELLSetPreallocation(Mat,PetscInt,const PetscInt[]);
322: PETSC_EXTERN PetscErrorCode MatMPISELLSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);

324: PETSC_EXTERN PetscErrorCode MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
325: PETSC_EXTERN PetscErrorCode MatCreateDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
326: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
327: PETSC_EXTERN PetscErrorCode MatCreateAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
328: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat*);
329: PETSC_EXTERN PetscErrorCode MatUpdateMPIAIJWithArrays(Mat,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
330: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*);
331: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithSeqAIJ(MPI_Comm,Mat,Mat,const PetscInt[],Mat*);

333: PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
334: PETSC_EXTERN PetscErrorCode MatCreateBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
335: PETSC_EXTERN PetscErrorCode MatCreateMPIBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat*);

337: PETSC_EXTERN PetscErrorCode MatSetPreallocationCOO(Mat,PetscCount,const PetscInt[],const PetscInt[]);
338: PETSC_EXTERN PetscErrorCode MatSetPreallocationCOOLocal(Mat,PetscCount,PetscInt[],PetscInt[]);
339: PETSC_EXTERN PetscErrorCode MatSetValuesCOO(Mat,const PetscScalar[],InsertMode);

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

344: PETSC_EXTERN PetscErrorCode MatCreateSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
345: PETSC_EXTERN PetscErrorCode MatCreateMPISBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
346: PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
347: PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
348: PETSC_EXTERN PetscErrorCode MatXAIJSetPreallocation(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscInt[],const PetscInt[]);

350: PETSC_EXTERN PetscErrorCode MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*);
351: PETSC_EXTERN PetscErrorCode MatCreateCentering(MPI_Comm,PetscInt,PetscInt,Mat*);
352: PETSC_EXTERN PetscErrorCode MatCreateNormal(Mat,Mat*);
353: PETSC_EXTERN PetscErrorCode MatCreateNormalHermitian(Mat,Mat*);
354: PETSC_EXTERN PetscErrorCode MatCreateLRC(Mat,Mat,Vec,Mat,Mat*);
355: PETSC_EXTERN PetscErrorCode MatLRCGetMats(Mat,Mat*,Mat*,Vec*,Mat*);
356: PETSC_EXTERN PetscErrorCode MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,ISLocalToGlobalMapping,Mat*);
357: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
358: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);

360: PETSC_EXTERN PetscErrorCode MatCreateScatter(MPI_Comm,VecScatter,Mat*);
361: PETSC_EXTERN PetscErrorCode MatScatterSetVecScatter(Mat,VecScatter);
362: PETSC_EXTERN PetscErrorCode MatScatterGetVecScatter(Mat,VecScatter*);
363: PETSC_EXTERN PetscErrorCode MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*);
364: PETSC_EXTERN PetscErrorCode MatCompositeAddMat(Mat,Mat);
365: PETSC_EXTERN PetscErrorCode MatCompositeMerge(Mat);
366: typedef enum {MAT_COMPOSITE_MERGE_RIGHT,MAT_COMPOSITE_MERGE_LEFT} MatCompositeMergeType;
367: PETSC_EXTERN PetscErrorCode MatCompositeSetMergeType(Mat,MatCompositeMergeType);
368: PETSC_EXTERN PetscErrorCode MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*);
369: typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType;
370: PETSC_EXTERN PetscErrorCode MatCompositeSetType(Mat,MatCompositeType);
371: PETSC_EXTERN PetscErrorCode MatCompositeGetType(Mat,MatCompositeType*);
372: PETSC_EXTERN PetscErrorCode MatCompositeSetMatStructure(Mat,MatStructure);
373: PETSC_EXTERN PetscErrorCode MatCompositeGetMatStructure(Mat,MatStructure*);
374: PETSC_EXTERN PetscErrorCode MatCompositeGetNumberMat(Mat,PetscInt*);
375: PETSC_EXTERN PetscErrorCode MatCompositeGetMat(Mat,PetscInt,Mat*);
376: PETSC_EXTERN PetscErrorCode MatCompositeSetScalings(Mat,const PetscScalar*);

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

381: PETSC_EXTERN PetscErrorCode MatCreateTranspose(Mat,Mat*);
382: PETSC_EXTERN PetscErrorCode MatTransposeGetMat(Mat,Mat*);
383: PETSC_EXTERN PetscErrorCode MatCreateHermitianTranspose(Mat,Mat*);
384: PETSC_EXTERN PetscErrorCode MatHermitianTransposeGetMat(Mat,Mat*);
385: PETSC_EXTERN PetscErrorCode MatNormalGetMat(Mat,Mat*);
386: PETSC_EXTERN PetscErrorCode MatNormalHermitianGetMat(Mat,Mat*);
387: PETSC_EXTERN PetscErrorCode MatCreateSubMatrixVirtual(Mat,IS,IS,Mat*);
388: PETSC_EXTERN PetscErrorCode MatSubMatrixVirtualUpdate(Mat,Mat,IS,IS);
389: PETSC_EXTERN PetscErrorCode MatCreateLocalRef(Mat,IS,IS,Mat*);
390: PETSC_EXTERN PetscErrorCode MatCreateConstantDiagonal(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar,Mat*);

392: #if defined(PETSC_HAVE_HYPRE)
393: PETSC_EXTERN PetscErrorCode MatHYPRESetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
394: #endif

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

398: PETSC_EXTERN PetscErrorCode MatResetPreallocation(Mat);
399: PETSC_EXTERN PetscErrorCode MatSetUp(Mat);
400: PETSC_EXTERN PetscErrorCode MatDestroy(Mat*);
401: PETSC_EXTERN PetscErrorCode MatGetNonzeroState(Mat,PetscObjectState*);

403: PETSC_EXTERN PetscErrorCode MatConjugate(Mat);
404: PETSC_EXTERN PetscErrorCode MatRealPart(Mat);
405: PETSC_EXTERN PetscErrorCode MatImaginaryPart(Mat);
406: PETSC_EXTERN PetscErrorCode MatGetDiagonalBlock(Mat,Mat*);
407: PETSC_EXTERN PetscErrorCode MatGetTrace(Mat,PetscScalar*);
408: PETSC_EXTERN PetscErrorCode MatInvertBlockDiagonal(Mat,const PetscScalar **);
409: PETSC_EXTERN PetscErrorCode MatInvertVariableBlockDiagonal(Mat,PetscInt,const PetscInt*,PetscScalar*);
410: PETSC_EXTERN PetscErrorCode MatInvertBlockDiagonalMat(Mat,Mat);

412: /* ------------------------------------------------------------*/
413: PETSC_EXTERN PetscErrorCode MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
414: PETSC_EXTERN PetscErrorCode MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
415: PETSC_EXTERN PetscErrorCode MatSetValuesRow(Mat,PetscInt,const PetscScalar[]);
416: PETSC_EXTERN PetscErrorCode MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]);
417: PETSC_EXTERN PetscErrorCode MatSetValuesBatch(Mat,PetscInt,PetscInt,PetscInt[],const PetscScalar[]);
418: PETSC_EXTERN PetscErrorCode MatSetRandom(Mat,PetscRandom);

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

424:    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).
425:    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.

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

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

431:    Level: beginner

433: .seealso:  MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockedStencil(), DMDAVecGetArray(), DMDAVecGetArrayF90()
434: S*/
435: typedef struct {
436:   PetscInt k,j,i,c;
437: } MatStencil;

439: PETSC_EXTERN PetscErrorCode MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
440: PETSC_EXTERN PetscErrorCode MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
441: PETSC_EXTERN PetscErrorCode MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt);

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

447:     Level: beginner

449: .seealso: MatAssemblyBegin(), MatAssemblyEnd()
450: E*/
451: typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType;
452: PETSC_EXTERN PetscErrorCode MatAssemblyBegin(Mat,MatAssemblyType);
453: PETSC_EXTERN PetscErrorCode MatAssemblyEnd(Mat,MatAssemblyType);
454: PETSC_EXTERN PetscErrorCode MatAssembled(Mat,PetscBool *);

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

459:     Level: beginner

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

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

467: .seealso: MatSetOption()
468: E*/
469: typedef enum {MAT_OPTION_MIN = -3,
470:               MAT_UNUSED_NONZERO_LOCATION_ERR = -2,
471:               MAT_ROW_ORIENTED = -1,
472:               MAT_SYMMETRIC = 1,
473:               MAT_STRUCTURALLY_SYMMETRIC = 2,
474:               MAT_FORCE_DIAGONAL_ENTRIES = 3,
475:               MAT_IGNORE_OFF_PROC_ENTRIES = 4,
476:               MAT_USE_HASH_TABLE = 5,
477:               MAT_KEEP_NONZERO_PATTERN = 6,
478:               MAT_IGNORE_ZERO_ENTRIES = 7,
479:               MAT_USE_INODES = 8,
480:               MAT_HERMITIAN = 9,
481:               MAT_SYMMETRY_ETERNAL = 10,
482:               MAT_NEW_NONZERO_LOCATION_ERR = 11,
483:               MAT_IGNORE_LOWER_TRIANGULAR = 12,
484:               MAT_ERROR_LOWER_TRIANGULAR = 13,
485:               MAT_GETROW_UPPERTRIANGULAR = 14,
486:               MAT_SPD = 15,
487:               MAT_NO_OFF_PROC_ZERO_ROWS = 16,
488:               MAT_NO_OFF_PROC_ENTRIES = 17,
489:               MAT_NEW_NONZERO_LOCATIONS = 18,
490:               MAT_NEW_NONZERO_ALLOCATION_ERR = 19,
491:               MAT_SUBSET_OFF_PROC_ENTRIES = 20,
492:               MAT_SUBMAT_SINGLEIS = 21,
493:               MAT_STRUCTURE_ONLY = 22,
494:               MAT_SORTED_FULL = 23,
495:               MAT_FORM_EXPLICIT_TRANSPOSE = 24,
496:               MAT_OPTION_MAX = 25} MatOption;

498: PETSC_EXTERN const char *const *MatOptions;
499: PETSC_EXTERN PetscErrorCode MatSetOption(Mat,MatOption,PetscBool);
500: PETSC_EXTERN PetscErrorCode MatGetOption(Mat,MatOption,PetscBool*);
501: PETSC_EXTERN PetscErrorCode MatPropagateSymmetryOptions(Mat,Mat);
502: PETSC_EXTERN PetscErrorCode MatGetType(Mat,MatType*);

504: PETSC_EXTERN PetscErrorCode MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
505: PETSC_EXTERN PetscErrorCode MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
506: PETSC_EXTERN PetscErrorCode MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
507: PETSC_EXTERN PetscErrorCode MatGetRowUpperTriangular(Mat);
508: PETSC_EXTERN PetscErrorCode MatRestoreRowUpperTriangular(Mat);
509: PETSC_EXTERN PetscErrorCode MatGetColumnVector(Mat,Vec,PetscInt);
510: PETSC_EXTERN PetscErrorCode MatSeqAIJGetArray(Mat,PetscScalar *[]);
511: PETSC_EXTERN PetscErrorCode MatSeqAIJGetArrayRead(Mat,const PetscScalar *[]);
512: PETSC_EXTERN PetscErrorCode MatSeqAIJGetArrayWrite(Mat,PetscScalar *[]);
513: PETSC_EXTERN PetscErrorCode MatSeqAIJRestoreArray(Mat,PetscScalar *[]);
514: PETSC_EXTERN PetscErrorCode MatSeqAIJRestoreArrayRead(Mat,const PetscScalar *[]);
515: PETSC_EXTERN PetscErrorCode MatSeqAIJRestoreArrayWrite(Mat,PetscScalar *[]);
516: PETSC_EXTERN PetscErrorCode MatSeqAIJGetMaxRowNonzeros(Mat,PetscInt*);
517: PETSC_EXTERN PetscErrorCode MatSeqAIJSetValuesLocalFast(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
518: PETSC_EXTERN PetscErrorCode MatSeqAIJSetType(Mat,MatType);
519: PETSC_EXTERN PetscErrorCode MatSeqAIJKron(Mat,Mat,MatReuse,Mat*);
520: PETSC_EXTERN PetscErrorCode MatSeqAIJRegister(const char[],PetscErrorCode (*)(Mat,MatType,MatReuse,Mat *));
521: PETSC_EXTERN PetscFunctionList MatSeqAIJList;
522: PETSC_EXTERN PetscErrorCode MatSeqBAIJGetArray(Mat,PetscScalar *[]);
523: PETSC_EXTERN PetscErrorCode MatSeqBAIJRestoreArray(Mat,PetscScalar *[]);
524: PETSC_EXTERN PetscErrorCode MatSeqSBAIJGetArray(Mat,PetscScalar *[]);
525: PETSC_EXTERN PetscErrorCode MatSeqSBAIJRestoreArray(Mat,PetscScalar *[]);
526: PETSC_EXTERN PetscErrorCode MatDenseGetArray(Mat,PetscScalar *[]);
527: PETSC_EXTERN PetscErrorCode MatDenseRestoreArray(Mat,PetscScalar *[]);
528: PETSC_EXTERN PetscErrorCode MatDensePlaceArray(Mat,const PetscScalar[]);
529: PETSC_EXTERN PetscErrorCode MatDenseReplaceArray(Mat,const PetscScalar[]);
530: PETSC_EXTERN PetscErrorCode MatDenseResetArray(Mat);
531: PETSC_EXTERN PetscErrorCode MatDenseGetArrayRead(Mat,const PetscScalar *[]);
532: PETSC_EXTERN PetscErrorCode MatDenseRestoreArrayRead(Mat,const PetscScalar *[]);
533: PETSC_EXTERN PetscErrorCode MatDenseGetArrayWrite(Mat,PetscScalar *[]);
534: PETSC_EXTERN PetscErrorCode MatDenseRestoreArrayWrite(Mat,PetscScalar *[]);
535: PETSC_EXTERN PetscErrorCode MatGetBlockSize(Mat,PetscInt *);
536: PETSC_EXTERN PetscErrorCode MatSetBlockSize(Mat,PetscInt);
537: PETSC_EXTERN PetscErrorCode MatGetBlockSizes(Mat,PetscInt *,PetscInt *);
538: PETSC_EXTERN PetscErrorCode MatSetBlockSizes(Mat,PetscInt,PetscInt);
539: PETSC_EXTERN PetscErrorCode MatSetBlockSizesFromMats(Mat,Mat,Mat);
540: PETSC_EXTERN PetscErrorCode MatSetVariableBlockSizes(Mat,PetscInt,PetscInt*);
541: PETSC_EXTERN PetscErrorCode MatGetVariableBlockSizes(Mat,PetscInt*,const PetscInt**);

543: PETSC_EXTERN PetscErrorCode MatDenseGetColumn(Mat,PetscInt,PetscScalar*[]);
544: PETSC_EXTERN PetscErrorCode MatDenseRestoreColumn(Mat,PetscScalar*[]);
545: PETSC_EXTERN PetscErrorCode MatDenseGetColumnVec(Mat,PetscInt,Vec*);
546: PETSC_EXTERN PetscErrorCode MatDenseRestoreColumnVec(Mat,PetscInt,Vec*);
547: PETSC_EXTERN PetscErrorCode MatDenseGetColumnVecRead(Mat,PetscInt,Vec*);
548: PETSC_EXTERN PetscErrorCode MatDenseRestoreColumnVecRead(Mat,PetscInt,Vec*);
549: PETSC_EXTERN PetscErrorCode MatDenseGetColumnVecWrite(Mat,PetscInt,Vec*);
550: PETSC_EXTERN PetscErrorCode MatDenseRestoreColumnVecWrite(Mat,PetscInt,Vec*);
551: PETSC_EXTERN PetscErrorCode MatDenseGetSubMatrix(Mat,PetscInt,PetscInt,Mat*);
552: PETSC_EXTERN PetscErrorCode MatDenseRestoreSubMatrix(Mat,Mat*);

554: PETSC_EXTERN PetscErrorCode MatMult(Mat,Vec,Vec);
555: PETSC_EXTERN PetscErrorCode MatMultDiagonalBlock(Mat,Vec,Vec);
556: PETSC_EXTERN PetscErrorCode MatMultAdd(Mat,Vec,Vec,Vec);
557: PETSC_EXTERN PetscErrorCode MatMultTranspose(Mat,Vec,Vec);
558: PETSC_EXTERN PetscErrorCode MatMultHermitianTranspose(Mat,Vec,Vec);
559: PETSC_EXTERN PetscErrorCode MatIsTranspose(Mat,Mat,PetscReal,PetscBool *);
560: PETSC_EXTERN PetscErrorCode MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *);
561: PETSC_EXTERN PetscErrorCode MatMultTransposeAdd(Mat,Vec,Vec,Vec);
562: PETSC_EXTERN PetscErrorCode MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec);
563: PETSC_EXTERN PetscErrorCode MatMatSolve(Mat,Mat,Mat);
564: PETSC_EXTERN PetscErrorCode MatMatSolveTranspose(Mat,Mat,Mat);
565: PETSC_EXTERN PetscErrorCode MatMatTransposeSolve(Mat,Mat,Mat);
566: PETSC_EXTERN PetscErrorCode MatResidual(Mat,Vec,Vec,Vec);

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

572:     Level: beginner

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

576: $   MAT_DO_NOT_COPY_VALUES    - Create a matrix using the same nonzero pattern as the original matrix,
577: $                               with zeros for the numerical values.
578: $   MAT_COPY_VALUES           - Create a matrix with the same nonzero pattern as the original matrix
579: $                               and with the same numerical values.
580: $   MAT_SHARE_NONZERO_PATTERN - Create a matrix that shares the nonzero structure with the previous matrix
581: $                               and does not copy it, using zeros for the numerical values. The parent and
582: $                               child matrices will share their index (i and j) arrays, and you cannot
583: $                               insert new nonzero entries into either matrix.

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

589: .seealso: MatDuplicate()
590: E*/
591: typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption;

593: PETSC_EXTERN PetscErrorCode MatConvert(Mat,MatType,MatReuse,Mat*);
594: PETSC_EXTERN PetscErrorCode MatDuplicate(Mat,MatDuplicateOption,Mat*);

596: PETSC_EXTERN PetscErrorCode MatCopy(Mat,Mat,MatStructure);
597: PETSC_EXTERN PetscErrorCode MatView(Mat,PetscViewer);
598: PETSC_EXTERN PetscErrorCode MatIsSymmetric(Mat,PetscReal,PetscBool *);
599: PETSC_EXTERN PetscErrorCode MatIsStructurallySymmetric(Mat,PetscBool *);
600: PETSC_EXTERN PetscErrorCode MatIsHermitian(Mat,PetscReal,PetscBool *);
601: PETSC_EXTERN PetscErrorCode MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *);
602: PETSC_EXTERN PetscErrorCode MatIsHermitianKnown(Mat,PetscBool *,PetscBool *);
603: PETSC_EXTERN PetscErrorCode MatMissingDiagonal(Mat,PetscBool  *,PetscInt *);
604: PETSC_EXTERN PetscErrorCode MatLoad(Mat, PetscViewer);

606: PETSC_EXTERN PetscErrorCode MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
607: PETSC_EXTERN PetscErrorCode MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool  *);
608: PETSC_EXTERN PetscErrorCode MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
609: PETSC_EXTERN PetscErrorCode MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool  *);

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

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

616:    Level: intermediate

618: .seealso:  MatGetInfo(), MatInfoType
619: S*/
620: typedef struct {
621:   PetscLogDouble block_size;                         /* block size */
622:   PetscLogDouble nz_allocated,nz_used,nz_unneeded;   /* number of nonzeros */
623:   PetscLogDouble memory;                             /* memory allocated */
624:   PetscLogDouble assemblies;                         /* number of matrix assemblies called */
625:   PetscLogDouble mallocs;                            /* number of mallocs during MatSetValues() */
626:   PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
627:   PetscLogDouble factor_mallocs;                     /* number of mallocs during factorization */
628: } MatInfo;

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

634:     Level: beginner

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

638: .seealso: MatGetInfo(), MatInfo
639: E*/
640: typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
641: PETSC_EXTERN PetscErrorCode MatGetInfo(Mat,MatInfoType,MatInfo*);
642: PETSC_EXTERN PetscErrorCode MatGetDiagonal(Mat,Vec);
643: PETSC_EXTERN PetscErrorCode MatGetRowMax(Mat,Vec,PetscInt[]);
644: PETSC_EXTERN PetscErrorCode MatGetRowMin(Mat,Vec,PetscInt[]);
645: PETSC_EXTERN PetscErrorCode MatGetRowMaxAbs(Mat,Vec,PetscInt[]);
646: PETSC_EXTERN PetscErrorCode MatGetRowMinAbs(Mat,Vec,PetscInt[]);
647: PETSC_EXTERN PetscErrorCode MatGetRowSum(Mat,Vec);
648: PETSC_EXTERN PetscErrorCode MatTranspose(Mat,MatReuse,Mat*);
649: PETSC_EXTERN PetscErrorCode MatHermitianTranspose(Mat,MatReuse,Mat*);
650: PETSC_EXTERN PetscErrorCode MatPermute(Mat,IS,IS,Mat*);
651: PETSC_EXTERN PetscErrorCode MatDiagonalScale(Mat,Vec,Vec);
652: PETSC_EXTERN PetscErrorCode MatDiagonalSet(Mat,Vec,InsertMode);

654: PETSC_EXTERN PetscErrorCode MatEqual(Mat,Mat,PetscBool*);
655: PETSC_EXTERN PetscErrorCode MatMultEqual(Mat,Mat,PetscInt,PetscBool*);
656: PETSC_EXTERN PetscErrorCode MatMultAddEqual(Mat,Mat,PetscInt,PetscBool*);
657: PETSC_EXTERN PetscErrorCode MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool*);
658: PETSC_EXTERN PetscErrorCode MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool*);
659: PETSC_EXTERN PetscErrorCode MatMultHermitianTransposeEqual(Mat,Mat,PetscInt,PetscBool*);
660: PETSC_EXTERN PetscErrorCode MatMultHermitianTransposeAddEqual(Mat,Mat,PetscInt,PetscBool*);
661: PETSC_EXTERN PetscErrorCode MatMatMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*);
662: PETSC_EXTERN PetscErrorCode MatTransposeMatMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*);
663: PETSC_EXTERN PetscErrorCode MatMatTransposeMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*);
664: PETSC_EXTERN PetscErrorCode MatPtAPMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*);
665: PETSC_EXTERN PetscErrorCode MatRARtMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*);
666: PETSC_EXTERN PetscErrorCode MatIsLinear(Mat,PetscInt,PetscBool*);

668: PETSC_EXTERN PetscErrorCode MatNorm(Mat,NormType,PetscReal*);
669: PETSC_EXTERN PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal*);
670: PETSC_EXTERN PetscErrorCode MatGetColumnSums(Mat,PetscScalar*);
671: PETSC_EXTERN PetscErrorCode MatGetColumnSumsRealPart(Mat,PetscReal*);
672: PETSC_EXTERN PetscErrorCode MatGetColumnSumsImaginaryPart(Mat,PetscReal*);
673: PETSC_EXTERN PetscErrorCode MatGetColumnMeans(Mat,PetscScalar*);
674: PETSC_EXTERN PetscErrorCode MatGetColumnMeansRealPart(Mat,PetscReal*);
675: PETSC_EXTERN PetscErrorCode MatGetColumnMeansImaginaryPart(Mat,PetscReal*);
676: PETSC_EXTERN PetscErrorCode MatGetColumnReductions(Mat,PetscInt,PetscReal*);
677: PETSC_EXTERN PetscErrorCode MatZeroEntries(Mat);
678: PETSC_EXTERN PetscErrorCode MatSetInf(Mat);
679: PETSC_EXTERN PetscErrorCode MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
680: PETSC_EXTERN PetscErrorCode MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec);
681: PETSC_EXTERN PetscErrorCode MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec);
682: PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsStencil(Mat,PetscInt,const MatStencil[],PetscScalar,Vec,Vec);
683: PETSC_EXTERN PetscErrorCode MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
684: PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec);

686: PETSC_EXTERN PetscErrorCode MatGetSize(Mat,PetscInt*,PetscInt*);
687: PETSC_EXTERN PetscErrorCode MatGetLocalSize(Mat,PetscInt*,PetscInt*);
688: PETSC_EXTERN PetscErrorCode MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
689: PETSC_EXTERN PetscErrorCode MatGetOwnershipRanges(Mat,const PetscInt**);
690: PETSC_EXTERN PetscErrorCode MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*);
691: PETSC_EXTERN PetscErrorCode MatGetOwnershipRangesColumn(Mat,const PetscInt**);
692: PETSC_EXTERN PetscErrorCode MatGetOwnershipIS(Mat,IS*,IS*);

694: PETSC_EXTERN PetscErrorCode MatCreateSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
695: PETSC_DEPRECATED_FUNCTION("Use MatCreateSubMatrices() (since version 3.8)") 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);}
696: PETSC_EXTERN PetscErrorCode MatCreateSubMatricesMPI(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
697: PETSC_DEPRECATED_FUNCTION("Use MatCreateSubMatricesMPI() (since version 3.8)") 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);}
698: PETSC_EXTERN PetscErrorCode MatDestroyMatrices(PetscInt,Mat *[]);
699: PETSC_EXTERN PetscErrorCode MatDestroySubMatrices(PetscInt,Mat *[]);
700: PETSC_EXTERN PetscErrorCode MatCreateSubMatrix(Mat,IS,IS,MatReuse,Mat *);
701: PETSC_DEPRECATED_FUNCTION("Use MatCreateSubMatrix() (since version 3.8)") static inline PetscErrorCode MatGetSubMatrix(Mat mat,IS isrow,IS iscol,MatReuse cll,Mat *newmat) {return MatCreateSubMatrix(mat,isrow,iscol,cll,newmat);}
702: PETSC_EXTERN PetscErrorCode MatGetLocalSubMatrix(Mat,IS,IS,Mat*);
703: PETSC_EXTERN PetscErrorCode MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*);
704: PETSC_EXTERN PetscErrorCode MatGetSeqNonzeroStructure(Mat,Mat*);
705: PETSC_EXTERN PetscErrorCode MatDestroySeqNonzeroStructure(Mat*);

707: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJ(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*);
708: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*);
709: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJNumeric(Mat,Mat);
710: PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMat(Mat,MatReuse,Mat*);
711: PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*);
712: PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMatMerge(Mat,MatReuse,IS*,Mat*);
713: PETSC_EXTERN PetscErrorCode MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,Mat*);
714: PETSC_EXTERN PetscErrorCode MatGetGhosts(Mat, PetscInt *,const PetscInt *[]);

716: PETSC_EXTERN PetscErrorCode MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt);
717: PETSC_EXTERN PetscErrorCode MatIncreaseOverlapSplit(Mat mat,PetscInt n,IS is[],PetscInt ov);
718: PETSC_EXTERN PetscErrorCode MatMPIAIJSetUseScalableIncreaseOverlap(Mat,PetscBool);

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

722: PETSC_EXTERN PetscErrorCode MatMatMatMult(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
723: PETSC_EXTERN PetscErrorCode MatGalerkin(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);

725: PETSC_EXTERN PetscErrorCode MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*);
726: PETSC_EXTERN PetscErrorCode MatRARt(Mat,Mat,MatReuse,PetscReal,Mat*);

728: PETSC_EXTERN PetscErrorCode MatTransposeMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
729: PETSC_EXTERN PetscErrorCode MatMatTransposeMult(Mat,Mat,MatReuse,PetscReal,Mat*);

731: PETSC_EXTERN PetscErrorCode MatAXPY(Mat,PetscScalar,Mat,MatStructure);
732: PETSC_EXTERN PetscErrorCode MatAYPX(Mat,PetscScalar,Mat,MatStructure);

734: PETSC_EXTERN PetscErrorCode MatScale(Mat,PetscScalar);
735: PETSC_EXTERN PetscErrorCode MatShift(Mat,PetscScalar);

737: PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
738: PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*);
739: PETSC_EXTERN PetscErrorCode MatGetLayouts(Mat,PetscLayout*,PetscLayout*);
740: PETSC_EXTERN PetscErrorCode MatSetLayouts(Mat,PetscLayout,PetscLayout);
741: PETSC_EXTERN PetscErrorCode MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
742: PETSC_EXTERN PetscErrorCode MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
743: PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
744: PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
745: PETSC_EXTERN PetscErrorCode MatGetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
746: PETSC_EXTERN PetscErrorCode MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
747: PETSC_EXTERN PetscErrorCode MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);

749: PETSC_EXTERN PetscErrorCode MatStashSetInitialSize(Mat,PetscInt,PetscInt);
750: PETSC_EXTERN PetscErrorCode MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);

752: PETSC_EXTERN PetscErrorCode MatInterpolate(Mat,Vec,Vec);
753: PETSC_EXTERN PetscErrorCode MatInterpolateAdd(Mat,Vec,Vec,Vec);
754: PETSC_EXTERN PetscErrorCode MatRestrict(Mat,Vec,Vec);
755: PETSC_EXTERN PetscErrorCode MatMatInterpolate(Mat,Mat,Mat*);
756: PETSC_EXTERN PetscErrorCode MatMatInterpolateAdd(Mat,Mat,Mat,Mat*);
757: PETSC_EXTERN PetscErrorCode MatMatRestrict(Mat,Mat,Mat*);
758: PETSC_EXTERN PetscErrorCode MatCreateVecs(Mat,Vec*,Vec*);
759: PETSC_DEPRECATED_FUNCTION("Use MatCreateVecs() (since version 3.6)") static inline PetscErrorCode MatGetVecs(Mat mat,Vec *x,Vec *y) {return MatCreateVecs(mat,x,y);}
760: PETSC_EXTERN PetscErrorCode MatCreateRedundantMatrix(Mat,PetscInt,MPI_Comm,MatReuse,Mat*);
761: PETSC_EXTERN PetscErrorCode MatGetMultiProcBlock(Mat,MPI_Comm,MatReuse,Mat*);
762: PETSC_EXTERN PetscErrorCode MatFindZeroDiagonals(Mat,IS*);
763: PETSC_EXTERN PetscErrorCode MatFindOffBlockDiagonalEntries(Mat,IS*);
764: PETSC_EXTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);

766: /*MC
767:    MatSetValue - Set a single entry into a matrix.

769:    Not collective

771:    Synopsis:
772: #include <petscmat.h>
773:      PetscErrorCode MatSetValue(Mat m,PetscInt row,PetscInt col,PetscScalar value,InsertMode mode)

775:    Input Parameters:
776: +  m - the matrix
777: .  row - the row location of the entry
778: .  col - the column location of the entry
779: .  value - the value to insert
780: -  mode - either INSERT_VALUES or ADD_VALUES

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

785:    Level: beginner

787: .seealso: MatGetValue(), MatSetValues(), MatSetValueLocal(), MatSetValuesLocal()
788: M*/
789: static inline PetscErrorCode MatSetValue(Mat v,PetscInt i,PetscInt j,PetscScalar va,InsertMode mode) {return MatSetValues(v,1,&i,1,&j,&va,mode);}

791: /*@C
792:    MatGetValue - Gets a single value from a matrix

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

796:    Input Parameters:
797: +  mat - the matrix
798: .  row - the row location of the entry
799: -  col - the column location of the entry

801:    Output Parameter:
802: .  va - the value

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

807:    See notes for MatGetValues().

809:    Level: advanced

811: .seealso: MatSetValue(), MatGetValueLocal(), MatGetValues()
812: @*/
813: static inline PetscErrorCode MatGetValue(Mat mat,PetscInt row,PetscInt col,PetscScalar *va) {return MatGetValues(mat,1,&row,1,&col,va);}

815: /*MC
816:    MatSetValueLocal - Inserts or adds a single value into a matrix,
817:    using a local numbering of the nodes.

819:    Not Collective

821:    Input Parameters:
822: +  m - the matrix
823: .  row - the row location of the entry
824: .  col - the column location of the entry
825: .  value - the value to insert
826: -  mode - either INSERT_VALUES or ADD_VALUES

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

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

833:    Level: intermediate

835: .seealso: MatSetValue(), MatSetValuesLocal()
836: M*/
837: static inline PetscErrorCode MatSetValueLocal(Mat v,PetscInt i,PetscInt j,PetscScalar va,InsertMode mode) {return MatSetValuesLocal(v,1,&i,1,&j,&va,mode);}

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

843:    Synopsis:
844: #include <petscmat.h>
845:    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)

847:    Collective

849:    Input Parameters:
850: +  comm - the communicator that will share the eventually allocated matrix
851: .  nrows - the number of LOCAL rows in the matrix
852: -  ncols - the number of LOCAL columns in the matrix

854:    Output Parameters:
855: +  dnz - the array that will be passed to the matrix preallocation routines
856: -  onz - the other array passed to the matrix preallocation routines

858:    Level: intermediate

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

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

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

867: .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
868:           MatPreallocateSymmetricSetLocalBlock()
869: M*/
870: #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; do {                             \
871:   PetscInt __nrows = (nrows),__ncols = (ncols),__rstart,__start,__end = 0;                     \
872:   PetscCalloc2(__nrows,&(dnz),__nrows,&(onz));                                        \
873:   MPI_Scan(&__ncols,&__end,1,MPIU_INT,MPI_SUM,comm);                                \
874:   __start = __end - __ncols; (void)__start;                                                    \
875:   MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);                             \
876:   __rstart -= __nrows

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

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

886:    Not Collective

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

898:    Level: intermediate

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

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

905: .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock()
906:           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocalRemoveDups()
907: M*/
908: #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz)                        \
909:   PetscMacroReturnStandard(                                                                    \
910:     ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);                                \
911:     ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);                                \
912:     for (PetscInt __l=0;__l<nrows;__l++) MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz); \
913:   )

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

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

923:    Not Collective

925:    Input Parameters:
926: +  map - the row mapping from local numbering to global numbering
927: .  nrows - the number of rows indicated
928: .  rows - the indices of the rows (these values are mapped to the global values)
929: .  cmap - the column mapping from local to global numbering
930: .  ncols - the number of columns in the matrix   (this value will be changed if duplicate columns are found)
931: .  cols - the columns indicated (these values are mapped to the global values, they are then sorted and duplicates removed)
932: .  dnz - the array that will be passed to the matrix preallocation routines
933: -  onz - the other array passed to the matrix preallocation routines

935:    Level: intermediate

937:    Notes:
938:     See Users-Manual: ch_performance for more details.

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

942: .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock()
943:           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal()
944: M*/
945: #define MatPreallocateSetLocalRemoveDups(rmap,nrows,rows,cmap,ncols,cols,dnz,onz)              \
946:   PetscMacroReturnStandard(                                                                    \
947:     ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);                                \
948:     ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);                                \
949:     PetscSortRemoveDupsInt(&ncols,cols);                                              \
950:     for (PetscInt __l=0;__l<nrows;__l++) MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz); \
951:   )

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

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

961:    Not Collective

963:    Input Parameters:
964: +  map - the row mapping from local numbering to global numbering
965: .  nrows - the number of rows indicated
966: .  rows - the indices of the rows
967: .  cmap - the column mapping from local to global numbering
968: .  ncols - the number of columns in the matrix
969: .  cols - the columns indicated
970: .  dnz - the array that will be passed to the matrix preallocation routines
971: -  onz - the other array passed to the matrix preallocation routines

973:    Level: intermediate

975:    Notes:
976:     See Users-Manual: ch_performance for more details.

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

980: .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock()
981:           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock()
982: M*/
983: #define MatPreallocateSetLocalBlock(rmap,nrows,rows,cmap,ncols,cols,dnz,onz)                   \
984:   PetscMacroReturnStandard(                                                                    \
985:     ISLocalToGlobalMappingApplyBlock(rmap,nrows,rows,rows);                           \
986:     ISLocalToGlobalMappingApplyBlock(cmap,ncols,cols,cols);                           \
987:     for (PetscInt __l=0;__l<nrows;__l++) MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz); \
988:   )

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

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

998:    Not Collective

1000:    Input Parameters:
1001: +  map - the mapping between local numbering and global numbering
1002: .  nrows - the number of rows indicated
1003: .  rows - the indices of the rows
1004: .  ncols - the number of columns in the matrix
1005: .  cols - the columns indicated
1006: .  dnz - the array that will be passed to the matrix preallocation routines
1007: -  onz - the other array passed to the matrix preallocation routines

1009:    Level: intermediate

1011:    Notes:
1012:     See Users-Manual: ch_performance for more details.

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

1016: .seealso: MatPreallocateFinalize(), MatPreallocateSet()
1017:           MatPreallocateInitialize(),  MatPreallocateSetLocal()
1018: M*/
1019: #define MatPreallocateSymmetricSetLocalBlock(map,nrows,rows,ncols,cols,dnz,onz)                \
1020:   PetscMacroReturnStandard(                                                                    \
1021:     ISLocalToGlobalMappingApplyBlock(map,nrows,rows,rows);                            \
1022:     ISLocalToGlobalMappingApplyBlock(map,ncols,cols,cols);                            \
1023:     for (PetscInt __l=0;__l<nrows;__l++) MatPreallocateSymmetricSetBlock((rows)[__l],ncols,cols,dnz,onz); \
1024:   )

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

1030:    Synopsis:
1031: #include <petscmat.h>
1032:    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)

1034:    Not Collective

1036:    Input Parameters:
1037: +  row - the row
1038: .  ncols - the number of columns in the matrix
1039: -  cols - the columns indicated

1041:    Output Parameters:
1042: +  dnz - the array that will be passed to the matrix preallocation routines
1043: -  onz - the other array passed to the matrix preallocation routines

1045:    Level: intermediate

1047:    Notes:
1048:     See Users-Manual: ch_performance for more details.

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

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

1054: .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock()
1055:           MatPreallocateInitialize(), MatPreallocateSetLocal()
1056: M*/
1057: #define MatPreallocateSet(row,nc,cols,dnz,onz) PetscMacroReturnStandard(                       \
1060:     for (PetscInt __i = 0; __i < nc; ++__i) {                                                  \
1061:       if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++;                \
1062:       else if (dnz[row - __rstart] < __ncols) dnz[row - __rstart]++;                           \
1063:     }                                                                                          \
1064:   )

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

1070:    Synopsis:
1071: #include <petscmat.h>
1072:    PetscErrorCode MatPreallocateSymmetricSetBlock(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)

1074:    Not Collective

1076:    Input Parameters:
1077: +  nrows - the number of rows indicated
1078: .  rows - the indices of the rows
1079: .  ncols - the number of columns in the matrix
1080: .  cols - the columns indicated
1081: .  dnz - the array that will be passed to the matrix preallocation routines
1082: -  onz - the other array passed to the matrix preallocation routines

1084:    Level: intermediate

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

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

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

1093: .seealso: MatPreallocateFinalize(), MatPreallocateSet(),  MatPreallocateInitialize(),
1094:           MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal()
1095: M*/
1096: #define MatPreallocateSymmetricSetBlock(row,nc,cols,dnz,onz)                                   \
1097:   PetscMacroReturnStandard(                                                                    \
1098:     for (PetscInt __i=0; __i<nc; __i++) {                                                      \
1099:       if (cols[__i] >= __end) onz[row - __rstart]++;                                           \
1100:       else if (cols[__i] >= row && dnz[row - __rstart] < __ncols) dnz[row - __rstart]++;       \
1101:     }                                                                                          \
1102:   )

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

1107:    Synopsis:
1108: #include <petscmat.h>
1109:    PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz)

1111:    Not Collective

1113:    Input Parameters:
1114: +  A - matrix
1115: .  row - row where values exist (must be local to this process)
1116: .  ncols - number of columns
1117: .  cols - columns with nonzeros
1118: .  dnz - the array that will be passed to the matrix preallocation routines
1119: -  onz - the other array passed to the matrix preallocation routines

1121:    Level: intermediate

1123:    Notes:
1124:     See Users-Manual: ch_performance for more details.

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

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

1130: .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
1131:           MatPreallocateSymmetricSetLocalBlock()
1132: M*/
1133: #define MatPreallocateLocation(A,row,ncols,cols,dnz,onz) PetscMacroReturnStandard(      \
1134:     if (A) MatSetValues(A,1,&row,ncols,cols,NULL,INSERT_VALUES);               \
1135:     else   MatPreallocateSet(row,ncols,cols,dnz,onz);                          \
1136:   )

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

1142:    Synopsis:
1143: #include <petscmat.h>
1144:    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)

1146:    Collective

1148:    Input Parameters:
1149: +  dnz - the array that was be passed to the matrix preallocation routines
1150: -  onz - the other array passed to the matrix preallocation routines

1152:    Level: intermediate

1154:    Notes:
1155:     See Users-Manual: ch_performance for more details.

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

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

1161: .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
1162:           MatPreallocateSymmetricSetLocalBlock()
1163: M*/
1164: #define MatPreallocateFinalize(dnz,onz) PetscFree2(dnz,onz);} while (0)

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

1169: PETSC_EXTERN PetscErrorCode MatInodeAdjustForInodes(Mat,IS*,IS*);
1170: PETSC_EXTERN PetscErrorCode MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);

1172: PETSC_EXTERN PetscErrorCode MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
1173: PETSC_EXTERN PetscErrorCode MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
1174: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
1175: PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
1176: PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
1177: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJFromTriple(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*,PetscInt,PetscBool);

1179: #define MAT_SKIP_ALLOCATION -4

1181: PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
1182: PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
1183: PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
1184: PETSC_EXTERN PetscErrorCode MatSeqAIJSetTotalPreallocation(Mat,PetscInt);

1186: PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1187: PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1188: PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1189: PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []);
1190: PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
1191: PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
1192: PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
1193: PETSC_EXTERN PetscErrorCode MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
1194: PETSC_EXTERN PetscErrorCode MatMPIAdjToSeq(Mat,Mat*);
1195: PETSC_EXTERN PetscErrorCode MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
1196: PETSC_EXTERN PetscErrorCode MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
1197: PETSC_EXTERN PetscErrorCode MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,const PetscInt*[]);
1198: PETSC_EXTERN PetscErrorCode MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,const PetscInt*[]);
1199: PETSC_EXTERN PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat,Mat*);

1201: PETSC_EXTERN PetscErrorCode MatDenseGetLDA(Mat,PetscInt*);
1202: PETSC_EXTERN PetscErrorCode MatDenseSetLDA(Mat,PetscInt);
1203: PETSC_DEPRECATED_FUNCTION("Use MatDenseSetLDA() (since version 3.14)") static inline PetscErrorCode MatSeqDenseSetLDA(Mat A,PetscInt lda) {return MatDenseSetLDA(A,lda);}
1204: PETSC_EXTERN PetscErrorCode MatDenseGetLocalMatrix(Mat,Mat*);

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

1208: PETSC_EXTERN PetscErrorCode MatStoreValues(Mat);
1209: PETSC_EXTERN PetscErrorCode MatRetrieveValues(Mat);

1211: PETSC_EXTERN PetscErrorCode MatFindNonzeroRows(Mat,IS*);
1212: PETSC_EXTERN PetscErrorCode MatFindZeroRows(Mat,IS*);
1213: /*
1214:   These routines are not usually accessed directly, rather solving is
1215:   done through the KSP and PC interfaces.
1216: */

1218: /*J
1219:     MatOrderingType - String with the name of a PETSc matrix ordering

1221:    Level: beginner

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

1226: .seealso: MatGetOrdering()
1227: J*/
1228: typedef const char* MatOrderingType;
1229: #define MATORDERINGNATURAL        "natural"
1230: #define MATORDERINGND             "nd"
1231: #define MATORDERING1WD            "1wd"
1232: #define MATORDERINGRCM            "rcm"
1233: #define MATORDERINGQMD            "qmd"
1234: #define MATORDERINGROWLENGTH      "rowlength"
1235: #define MATORDERINGWBM            "wbm"
1236: #define MATORDERINGSPECTRAL       "spectral"
1237: #define MATORDERINGAMD            "amd"            /* only works if UMFPACK is installed with PETSc */
1238: #define MATORDERINGMETISND        "metisnd"        /* only works if METIS is installed with PETSc */
1239: #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 */
1240: #define MATORDERINGEXTERNAL       "external"       /* uses an ordering type internal to the factorization package */

1242: PETSC_EXTERN PetscErrorCode MatGetOrdering(Mat,MatOrderingType,IS*,IS*);
1243: PETSC_EXTERN PetscErrorCode MatGetOrderingList(PetscFunctionList*);
1244: PETSC_EXTERN PetscErrorCode MatOrderingRegister(const char[],PetscErrorCode(*)(Mat,MatOrderingType,IS*,IS*));
1245: PETSC_EXTERN PetscFunctionList MatOrderingList;

1247: PETSC_EXTERN PetscErrorCode MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
1248: PETSC_EXTERN PetscErrorCode MatCreateLaplacian(Mat,PetscReal,PetscBool,Mat*);

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

1252: /*S
1253:     MatFactorShiftType - Numeric Shift for factorizations

1255:    Level: beginner

1257: .seealso: MatGetFactor()
1258: S*/
1259: typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType;
1260: PETSC_EXTERN const char *const MatFactorShiftTypes[];
1261: PETSC_EXTERN const char *const MatFactorShiftTypesDetail[];

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

1266:     Level: beginner

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

1271: .seealso: MatGetFactor()
1272: S*/
1273: typedef enum {MAT_FACTOR_NOERROR,MAT_FACTOR_STRUCT_ZEROPIVOT,MAT_FACTOR_NUMERIC_ZEROPIVOT,MAT_FACTOR_OUTMEMORY,MAT_FACTOR_OTHER} MatFactorError;

1275: PETSC_EXTERN PetscErrorCode MatFactorGetError(Mat,MatFactorError*);
1276: PETSC_EXTERN PetscErrorCode MatFactorClearError(Mat);
1277: PETSC_EXTERN PetscErrorCode MatFactorGetErrorZeroPivot(Mat,PetscReal*,PetscInt*);

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

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

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

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

1290:    Level: developer

1292: .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
1293:           MatFactorInfoInitialize()

1295: S*/
1296: typedef struct {
1297:   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
1298:   PetscReal     usedt;
1299:   PetscReal     dt;             /* drop tolerance */
1300:   PetscReal     dtcol;          /* tolerance for pivoting */
1301:   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
1302:   PetscReal     fill;           /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */
1303:   PetscReal     levels;         /* ICC/ILU(levels) */
1304:   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
1305:                                    factorization may be faster if do not pivot */
1306:   PetscReal     zeropivot;      /* pivot is called zero if less than this */
1307:   PetscReal     shifttype;      /* type of shift added to matrix factor to prevent zero pivots */
1308:   PetscReal     shiftamount;     /* how large the shift is */
1309: } MatFactorInfo;

1311: PETSC_EXTERN PetscErrorCode MatFactorInfoInitialize(MatFactorInfo*);
1312: PETSC_EXTERN PetscErrorCode MatCholeskyFactor(Mat,IS,const MatFactorInfo*);
1313: PETSC_EXTERN PetscErrorCode MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
1314: PETSC_EXTERN PetscErrorCode MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*);
1315: PETSC_EXTERN PetscErrorCode MatLUFactor(Mat,IS,IS,const MatFactorInfo*);
1316: PETSC_EXTERN PetscErrorCode MatILUFactor(Mat,IS,IS,const MatFactorInfo*);
1317: PETSC_EXTERN PetscErrorCode MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
1318: PETSC_EXTERN PetscErrorCode MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
1319: PETSC_EXTERN PetscErrorCode MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
1320: PETSC_EXTERN PetscErrorCode MatICCFactor(Mat,IS,const MatFactorInfo*);
1321: PETSC_EXTERN PetscErrorCode MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*);
1322: PETSC_EXTERN PetscErrorCode MatQRFactor(Mat,IS,const MatFactorInfo*);
1323: PETSC_EXTERN PetscErrorCode MatQRFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
1324: PETSC_EXTERN PetscErrorCode MatQRFactorNumeric(Mat,Mat,const MatFactorInfo*);
1325: PETSC_EXTERN PetscErrorCode MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
1326: PETSC_EXTERN PetscErrorCode MatSolve(Mat,Vec,Vec);
1327: PETSC_EXTERN PetscErrorCode MatForwardSolve(Mat,Vec,Vec);
1328: PETSC_EXTERN PetscErrorCode MatBackwardSolve(Mat,Vec,Vec);
1329: PETSC_EXTERN PetscErrorCode MatSolveAdd(Mat,Vec,Vec,Vec);
1330: PETSC_EXTERN PetscErrorCode MatSolveTranspose(Mat,Vec,Vec);
1331: PETSC_EXTERN PetscErrorCode MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
1332: PETSC_EXTERN PetscErrorCode MatSolves(Mat,Vecs,Vecs);
1333: PETSC_EXTERN PetscErrorCode MatSetUnfactored(Mat);

1335: typedef enum {MAT_FACTOR_SCHUR_UNFACTORED, MAT_FACTOR_SCHUR_FACTORED, MAT_FACTOR_SCHUR_INVERTED} MatFactorSchurStatus;
1336: PETSC_EXTERN PetscErrorCode MatFactorSetSchurIS(Mat,IS);
1337: PETSC_EXTERN PetscErrorCode MatFactorGetSchurComplement(Mat,Mat*,MatFactorSchurStatus*);
1338: PETSC_EXTERN PetscErrorCode MatFactorRestoreSchurComplement(Mat,Mat*,MatFactorSchurStatus);
1339: PETSC_EXTERN PetscErrorCode MatFactorInvertSchurComplement(Mat);
1340: PETSC_EXTERN PetscErrorCode MatFactorCreateSchurComplement(Mat,Mat*,MatFactorSchurStatus*);
1341: PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplement(Mat,Vec,Vec);
1342: PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplementTranspose(Mat,Vec,Vec);
1343: PETSC_EXTERN PetscErrorCode MatFactorFactorizeSchurComplement(Mat);

1345: /*E
1346:     MatSORType - What type of (S)SOR to perform

1348:     Level: beginner

1350:    May be bitwise ORd together

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

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

1356: .seealso: MatSOR()
1357: E*/
1358: typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1359:               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1360:               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
1361:               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
1362: PETSC_EXTERN PetscErrorCode MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);

1364: /*
1365:     These routines are for efficiently computing Jacobians via finite differences.
1366: */

1368: /*S
1369:      MatColoring - Object for managing the coloring of matrices.

1371:    Level: beginner

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

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

1380: .seealso:  MatFDColoringCreate(), MatColoringWeightType, ISColoring, MatFDColoring, DMCreateColoring(), MatColoringCreate(), MatOrdering, MatPartitioning, MatColoringType
1381: S*/
1382: typedef struct _p_MatColoring* MatColoring;

1384: /*J
1385:     MatColoringType - String with the name of a PETSc matrix coloring

1387:    Level: beginner

1389: .seealso: MatColoringSetType(), MatColoring
1390: J*/
1391: typedef const  char*           MatColoringType;
1392: #define MATCOLORINGJP      "jp"
1393: #define MATCOLORINGPOWER   "power"
1394: #define MATCOLORINGNATURAL "natural"
1395: #define MATCOLORINGSL      "sl"
1396: #define MATCOLORINGLF      "lf"
1397: #define MATCOLORINGID      "id"
1398: #define MATCOLORINGGREEDY  "greedy"

1400: /*E
1401:    MatColoringWeightType - Type of weight scheme

1403:     Not Collective

1405: +   MAT_COLORING_RANDOM  - Random weights
1406: .   MAT_COLORING_LEXICAL - Lexical weighting based upon global numbering.
1407: -   MAT_COLORING_LF      - Last-first weighting.

1409:     Level: intermediate

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

1413: .seealso: MatColoring, MatColoringCreate()
1414: E*/
1415: typedef enum {MAT_COLORING_WEIGHT_RANDOM,MAT_COLORING_WEIGHT_LEXICAL,MAT_COLORING_WEIGHT_LF,MAT_COLORING_WEIGHT_SL} MatColoringWeightType;

1417: PETSC_EXTERN PetscErrorCode MatColoringCreate(Mat,MatColoring*);
1418: PETSC_EXTERN PetscErrorCode MatColoringGetDegrees(Mat,PetscInt,PetscInt*);
1419: PETSC_EXTERN PetscErrorCode MatColoringDestroy(MatColoring*);
1420: PETSC_EXTERN PetscErrorCode MatColoringView(MatColoring,PetscViewer);
1421: PETSC_EXTERN PetscErrorCode MatColoringSetType(MatColoring,MatColoringType);
1422: PETSC_EXTERN PetscErrorCode MatColoringSetFromOptions(MatColoring);
1423: PETSC_EXTERN PetscErrorCode MatColoringSetDistance(MatColoring,PetscInt);
1424: PETSC_EXTERN PetscErrorCode MatColoringGetDistance(MatColoring,PetscInt*);
1425: PETSC_EXTERN PetscErrorCode MatColoringSetMaxColors(MatColoring,PetscInt);
1426: PETSC_EXTERN PetscErrorCode MatColoringGetMaxColors(MatColoring,PetscInt*);
1427: PETSC_EXTERN PetscErrorCode MatColoringApply(MatColoring,ISColoring*);
1428: PETSC_EXTERN PetscErrorCode MatColoringRegister(const char[],PetscErrorCode(*)(MatColoring));
1429: PETSC_EXTERN PetscErrorCode MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
1430: PETSC_EXTERN PetscErrorCode MatColoringSetWeightType(MatColoring,MatColoringWeightType);
1431: PETSC_EXTERN PetscErrorCode MatColoringSetWeights(MatColoring,PetscReal*,PetscInt*);
1432: PETSC_EXTERN PetscErrorCode MatColoringCreateWeights(MatColoring,PetscReal **,PetscInt **lperm);
1433: PETSC_EXTERN PetscErrorCode MatColoringTest(MatColoring,ISColoring);
1434: PETSC_DEPRECATED_FUNCTION("Use MatColoringTest() (since version 3.10)") static inline PetscErrorCode MatColoringTestValid(MatColoring matcoloring,ISColoring iscoloring) {return MatColoringTest(matcoloring,iscoloring);}
1435: PETSC_EXTERN PetscErrorCode MatISColoringTest(Mat,ISColoring);

1437: /*S
1438:      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1439:         and coloring

1441:    Level: beginner

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

1446: .seealso:  MatFDColoringCreate(), MatColoring, DMCreateColoring()
1447: S*/
1448: typedef struct _p_MatFDColoring* MatFDColoring;

1450: PETSC_EXTERN PetscErrorCode MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
1451: PETSC_EXTERN PetscErrorCode MatFDColoringDestroy(MatFDColoring*);
1452: PETSC_EXTERN PetscErrorCode MatFDColoringView(MatFDColoring,PetscViewer);
1453: PETSC_EXTERN PetscErrorCode MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
1454: PETSC_EXTERN PetscErrorCode MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**);
1455: PETSC_EXTERN PetscErrorCode MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
1456: PETSC_EXTERN PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring);
1457: PETSC_EXTERN PetscErrorCode MatFDColoringApply(Mat,MatFDColoring,Vec,void *);
1458: PETSC_EXTERN PetscErrorCode MatFDColoringSetF(MatFDColoring,Vec);
1459: PETSC_EXTERN PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,const PetscInt*[]);
1460: PETSC_EXTERN PetscErrorCode MatFDColoringSetUp(Mat,ISColoring,MatFDColoring);
1461: PETSC_EXTERN PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring,PetscInt,PetscInt);
1462: PETSC_EXTERN PetscErrorCode MatFDColoringSetValues(Mat,MatFDColoring,const PetscScalar*);

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

1467:    Level: beginner

1469: .seealso:  MatTransposeColoringCreate()
1470: S*/
1471: typedef struct _p_MatTransposeColoring* MatTransposeColoring;

1473: PETSC_EXTERN PetscErrorCode MatTransposeColoringCreate(Mat,ISColoring,MatTransposeColoring *);
1474: PETSC_EXTERN PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring,Mat,Mat);
1475: PETSC_EXTERN PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring,Mat,Mat);
1476: PETSC_EXTERN PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring*);

1478: /*
1479:     These routines are for partitioning matrices: currently used only
1480:   for adjacency matrix, MatCreateMPIAdj().
1481: */

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

1486:    Level: beginner

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

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

1495: .seealso:  MatPartitioningCreate(), MatPartitioningType, MatColoring, MatGetOrdering()
1496: S*/
1497: typedef struct _p_MatPartitioning* MatPartitioning;

1499: /*J
1500:     MatPartitioningType - String with the name of a PETSc matrix partitioning

1502:    Level: beginner
1503: dm
1504: .seealso: MatPartitioningCreate(), MatPartitioning
1505: J*/
1506: typedef const char* MatPartitioningType;
1507: #define MATPARTITIONINGCURRENT  "current"
1508: #define MATPARTITIONINGAVERAGE   "average"
1509: #define MATPARTITIONINGSQUARE   "square"
1510: #define MATPARTITIONINGPARMETIS "parmetis"
1511: #define MATPARTITIONINGCHACO    "chaco"
1512: #define MATPARTITIONINGPARTY    "party"
1513: #define MATPARTITIONINGPTSCOTCH "ptscotch"
1514: #define MATPARTITIONINGHIERARCH  "hierarch"

1516: PETSC_EXTERN PetscErrorCode MatPartitioningCreate(MPI_Comm,MatPartitioning*);
1517: PETSC_EXTERN PetscErrorCode MatPartitioningSetType(MatPartitioning,MatPartitioningType);
1518: PETSC_EXTERN PetscErrorCode MatPartitioningSetNParts(MatPartitioning,PetscInt);
1519: PETSC_EXTERN PetscErrorCode MatPartitioningSetAdjacency(MatPartitioning,Mat);
1520: PETSC_EXTERN PetscErrorCode MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
1521: PETSC_EXTERN PetscErrorCode MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
1522: PETSC_EXTERN PetscErrorCode MatPartitioningSetUseEdgeWeights(MatPartitioning,PetscBool);
1523: PETSC_EXTERN PetscErrorCode MatPartitioningGetUseEdgeWeights(MatPartitioning,PetscBool*);
1524: PETSC_EXTERN PetscErrorCode MatPartitioningApply(MatPartitioning,IS*);
1525: PETSC_EXTERN PetscErrorCode MatPartitioningImprove(MatPartitioning,IS*);
1526: PETSC_EXTERN PetscErrorCode MatPartitioningViewImbalance(MatPartitioning,IS);
1527: PETSC_EXTERN PetscErrorCode MatPartitioningApplyND(MatPartitioning,IS*);
1528: PETSC_EXTERN PetscErrorCode MatPartitioningDestroy(MatPartitioning*);
1529: PETSC_EXTERN PetscErrorCode MatPartitioningRegister(const char[],PetscErrorCode (*)(MatPartitioning));
1530: PETSC_EXTERN PetscErrorCode MatPartitioningView(MatPartitioning,PetscViewer);
1531: PETSC_EXTERN PetscErrorCode MatPartitioningViewFromOptions(MatPartitioning,PetscObject,const char[]);
1532: PETSC_EXTERN PetscErrorCode MatPartitioningSetFromOptions(MatPartitioning);
1533: PETSC_EXTERN PetscErrorCode MatPartitioningGetType(MatPartitioning,MatPartitioningType*);

1535: PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetRepartition(MatPartitioning part);
1536: PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
1537: PETSC_EXTERN PetscErrorCode MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *);

1539: typedef enum { MP_CHACO_MULTILEVEL=1,MP_CHACO_SPECTRAL=2,MP_CHACO_LINEAR=4,MP_CHACO_RANDOM=5,MP_CHACO_SCATTERED=6 } MPChacoGlobalType;
1540: PETSC_EXTERN const char *const MPChacoGlobalTypes[];
1541: typedef enum { MP_CHACO_KERNIGHAN=1,MP_CHACO_NONE=2 } MPChacoLocalType;
1542: PETSC_EXTERN const char *const MPChacoLocalTypes[];
1543: typedef enum { MP_CHACO_LANCZOS=0,MP_CHACO_RQI=1 } MPChacoEigenType;
1544: PETSC_EXTERN const char *const MPChacoEigenTypes[];

1546: PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetGlobal(MatPartitioning,MPChacoGlobalType);
1547: PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetGlobal(MatPartitioning,MPChacoGlobalType*);
1548: PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetLocal(MatPartitioning,MPChacoLocalType);
1549: PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetLocal(MatPartitioning,MPChacoLocalType*);
1550: PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
1551: PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
1552: PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenSolver(MatPartitioning,MPChacoEigenType*);
1553: PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenTol(MatPartitioning,PetscReal);
1554: PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenTol(MatPartitioning,PetscReal*);
1555: PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenNumber(MatPartitioning,PetscInt);
1556: PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenNumber(MatPartitioning,PetscInt*);

1558: #define MP_PARTY_OPT "opt"
1559: #define MP_PARTY_LIN "lin"
1560: #define MP_PARTY_SCA "sca"
1561: #define MP_PARTY_RAN "ran"
1562: #define MP_PARTY_GBF "gbf"
1563: #define MP_PARTY_GCF "gcf"
1564: #define MP_PARTY_BUB "bub"
1565: #define MP_PARTY_DEF "def"
1566: PETSC_EXTERN PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning,const char*);
1567: #define MP_PARTY_HELPFUL_SETS "hs"
1568: #define MP_PARTY_KERNIGHAN_LIN "kl"
1569: #define MP_PARTY_NONE "no"
1570: PETSC_EXTERN PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning,const char*);
1571: PETSC_EXTERN PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
1572: PETSC_EXTERN PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning,PetscBool);
1573: PETSC_EXTERN PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool);

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

1578: PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning,PetscReal);
1579: PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning,PetscReal*);
1580: PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning,MPPTScotchStrategyType);
1581: PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning,MPPTScotchStrategyType*);

1583: /*
1584:  * hierarchical partitioning
1585:  */
1586: PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetFineparts(MatPartitioning,IS*);
1587: PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetCoarseparts(MatPartitioning,IS*);
1588: PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNcoarseparts(MatPartitioning,PetscInt);
1589: PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNfineparts(MatPartitioning, PetscInt);

1591: PETSC_EXTERN PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*);
1592: PETSC_EXTERN PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*);

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

1761: /*
1762:    Codes for matrices stored on disk. By default they are
1763:    stored in a universal format. By changing the format with
1764:    PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will
1765:    be stored in a way natural for the matrix, for example dense matrices
1766:    would be stored as dense. Matrices stored this way may only be
1767:    read into matrices of the same type.
1768: */
1769: #define MATRIX_BINARY_FORMAT_DENSE -1

1771: PETSC_EXTERN PetscErrorCode MatMPIBAIJSetHashTableFactor(Mat,PetscReal);

1773: PETSC_EXTERN PetscErrorCode MatISSetLocalMatType(Mat,MatType);
1774: PETSC_EXTERN PetscErrorCode MatISSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1775: PETSC_EXTERN PetscErrorCode MatISStoreL2L(Mat,PetscBool);
1776: PETSC_EXTERN PetscErrorCode MatISFixLocalEmpty(Mat,PetscBool);
1777: PETSC_EXTERN PetscErrorCode MatISGetLocalMat(Mat,Mat*);
1778: PETSC_EXTERN PetscErrorCode MatISRestoreLocalMat(Mat,Mat*);
1779: PETSC_EXTERN PetscErrorCode MatISSetLocalMat(Mat,Mat);
1780: PETSC_EXTERN PetscErrorCode MatISGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*);
1781: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use the MatConvert() interface (since version 3.10)") PetscErrorCode MatISGetMPIXAIJ(Mat,MatReuse,Mat*);

1783: /*S
1784:      MatNullSpace - Object that removes a null space from a vector, i.e.
1785:          orthogonalizes the vector to a subspace

1787:    Level: advanced

1789: .seealso:  MatNullSpaceCreate()
1790: S*/
1791: typedef struct _p_MatNullSpace* MatNullSpace;

1793: PETSC_EXTERN PetscErrorCode MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*);
1794: PETSC_EXTERN PetscErrorCode MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*);
1795: PETSC_EXTERN PetscErrorCode MatNullSpaceDestroy(MatNullSpace*);
1796: PETSC_EXTERN PetscErrorCode MatNullSpaceRemove(MatNullSpace,Vec);
1797: PETSC_EXTERN PetscErrorCode MatGetNullSpace(Mat, MatNullSpace *);
1798: PETSC_EXTERN PetscErrorCode MatGetTransposeNullSpace(Mat, MatNullSpace *);
1799: PETSC_EXTERN PetscErrorCode MatSetTransposeNullSpace(Mat,MatNullSpace);
1800: PETSC_EXTERN PetscErrorCode MatSetNullSpace(Mat,MatNullSpace);
1801: PETSC_EXTERN PetscErrorCode MatSetNearNullSpace(Mat,MatNullSpace);
1802: PETSC_EXTERN PetscErrorCode MatGetNearNullSpace(Mat,MatNullSpace*);
1803: PETSC_EXTERN PetscErrorCode MatNullSpaceTest(MatNullSpace,Mat,PetscBool  *);
1804: PETSC_EXTERN PetscErrorCode MatNullSpaceView(MatNullSpace,PetscViewer);
1805: PETSC_EXTERN PetscErrorCode MatNullSpaceGetVecs(MatNullSpace,PetscBool*,PetscInt*,const Vec**);
1806: PETSC_EXTERN PetscErrorCode MatNullSpaceCreateRigidBody(Vec,MatNullSpace*);

1808: PETSC_EXTERN PetscErrorCode MatReorderingSeqSBAIJ(Mat,IS);
1809: PETSC_EXTERN PetscErrorCode MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
1810: PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);

1812: PETSC_EXTERN PetscErrorCode MatCreateMAIJ(Mat,PetscInt,Mat*);
1813: PETSC_EXTERN PetscErrorCode MatMAIJRedimension(Mat,PetscInt,Mat*);
1814: PETSC_EXTERN PetscErrorCode MatMAIJGetAIJ(Mat,Mat*);

1816: PETSC_EXTERN PetscErrorCode MatComputeOperator(Mat,MatType,Mat*);
1817: PETSC_EXTERN PetscErrorCode MatComputeOperatorTranspose(Mat,MatType,Mat*);

1819: PETSC_DEPRECATED_FUNCTION("Use MatComputeOperator() (since version 3.12)") static inline PetscErrorCode MatComputeExplicitOperator(Mat A,Mat* B) { return MatComputeOperator(A,NULL,B); }
1820: PETSC_DEPRECATED_FUNCTION("Use MatComputeOperatorTranspose() (since version 3.12)") static inline PetscErrorCode MatComputeExplicitOperatorTranspose(Mat A,Mat* B) { return MatComputeOperatorTranspose(A,NULL,B); }

1822: PETSC_EXTERN PetscErrorCode MatCreateKAIJ(Mat,PetscInt,PetscInt,const PetscScalar[],const PetscScalar[],Mat*);
1823: PETSC_EXTERN PetscErrorCode MatKAIJGetAIJ(Mat,Mat*);
1824: PETSC_EXTERN PetscErrorCode MatKAIJGetS(Mat,PetscInt*,PetscInt*,PetscScalar**);
1825: PETSC_EXTERN PetscErrorCode MatKAIJGetSRead(Mat,PetscInt*,PetscInt*,const PetscScalar**);
1826: PETSC_EXTERN PetscErrorCode MatKAIJRestoreS(Mat,PetscScalar**);
1827: PETSC_EXTERN PetscErrorCode MatKAIJRestoreSRead(Mat,const PetscScalar**);
1828: PETSC_EXTERN PetscErrorCode MatKAIJGetT(Mat,PetscInt*,PetscInt*,PetscScalar**);
1829: PETSC_EXTERN PetscErrorCode MatKAIJGetTRead(Mat,PetscInt*,PetscInt*,const PetscScalar**);
1830: PETSC_EXTERN PetscErrorCode MatKAIJRestoreT(Mat,PetscScalar**);
1831: PETSC_EXTERN PetscErrorCode MatKAIJRestoreTRead(Mat,const PetscScalar**);
1832: PETSC_EXTERN PetscErrorCode MatKAIJSetAIJ(Mat,Mat);
1833: PETSC_EXTERN PetscErrorCode MatKAIJSetS(Mat,PetscInt,PetscInt,const PetscScalar[]);
1834: PETSC_EXTERN PetscErrorCode MatKAIJSetT(Mat,PetscInt,PetscInt,const PetscScalar[]);
1835: PETSC_EXTERN PetscErrorCode MatKAIJGetScaledIdentity(Mat,PetscBool*);

1837: PETSC_EXTERN PetscErrorCode MatDiagonalScaleLocal(Mat,Vec);

1839: PETSC_EXTERN PetscErrorCode MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
1840: PETSC_EXTERN PetscErrorCode MatMFFDSetBase(Mat,Vec,Vec);
1841: PETSC_EXTERN PetscErrorCode MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*);
1842: PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*));
1843: PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec));
1844: PETSC_EXTERN PetscErrorCode MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt);
1845: PETSC_EXTERN PetscErrorCode MatMFFDResetHHistory(Mat);
1846: PETSC_EXTERN PetscErrorCode MatMFFDSetFunctionError(Mat,PetscReal);
1847: PETSC_EXTERN PetscErrorCode MatMFFDSetPeriod(Mat,PetscInt);
1848: PETSC_EXTERN PetscErrorCode MatMFFDGetH(Mat,PetscScalar *);
1849: PETSC_EXTERN PetscErrorCode MatMFFDSetOptionsPrefix(Mat,const char[]);
1850: PETSC_EXTERN PetscErrorCode MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*);
1851: PETSC_EXTERN PetscErrorCode MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*);

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

1857:     Notes:
1858:     MATMFFD is a specific MatType which uses the MatMFFD data structure

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

1862:     Level: developer

1864: .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister()
1865: S*/
1866: typedef struct _p_MatMFFD* MatMFFD;

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

1871:    Level: beginner

1873: .seealso: MatMFFDSetType(), MatMFFDRegister()
1874: J*/
1875: typedef const char* MatMFFDType;
1876: #define MATMFFD_DS  "ds"
1877: #define MATMFFD_WP  "wp"

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

1882: PETSC_EXTERN PetscErrorCode MatMFFDDSSetUmin(Mat,PetscReal);
1883: PETSC_EXTERN PetscErrorCode MatMFFDWPSetComputeNormU(Mat,PetscBool);

1885: PETSC_EXTERN PetscErrorCode MatFDColoringSetType(MatFDColoring,MatMFFDType);

1887: PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
1888: PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);

1890: #ifdef PETSC_HAVE_H2OPUS
1891: PETSC_EXTERN_TYPEDEF typedef PetscScalar (*MatH2OpusKernel)(PetscInt,PetscReal[],PetscReal[],void*);
1892: PETSC_EXTERN PetscErrorCode MatCreateH2OpusFromKernel(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscReal[],PetscBool,MatH2OpusKernel,void*,PetscReal,PetscInt,PetscInt,Mat*);
1893: PETSC_EXTERN PetscErrorCode MatCreateH2OpusFromMat(Mat,PetscInt,const PetscReal[],PetscBool,PetscReal,PetscInt,PetscInt,PetscInt,PetscReal,Mat*);
1894: PETSC_EXTERN PetscErrorCode MatH2OpusSetSamplingMat(Mat,Mat,PetscInt,PetscReal);
1895: PETSC_EXTERN PetscErrorCode MatH2OpusOrthogonalize(Mat);
1896: PETSC_EXTERN PetscErrorCode MatH2OpusCompress(Mat,PetscReal);
1897: PETSC_EXTERN PetscErrorCode MatH2OpusSetNativeMult(Mat,PetscBool);
1898: PETSC_EXTERN PetscErrorCode MatH2OpusGetNativeMult(Mat,PetscBool*);
1899: PETSC_EXTERN PetscErrorCode MatH2OpusGetIndexMap(Mat,IS*);
1900: PETSC_EXTERN PetscErrorCode MatH2OpusMapVec(Mat,PetscBool,Vec,Vec*);
1901: PETSC_EXTERN PetscErrorCode MatH2OpusLowRankUpdate(Mat,Mat,Mat,PetscScalar);
1902: #endif

1904: #ifdef PETSC_HAVE_HTOOL
1905: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*MatHtoolKernel)(PetscInt,PetscInt,PetscInt,const PetscInt*,const PetscInt*,PetscScalar*,void*);
1906: PETSC_EXTERN PetscErrorCode MatCreateHtoolFromKernel(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscReal[],const PetscReal[],MatHtoolKernel,void*,Mat*);
1907: PETSC_EXTERN PetscErrorCode MatHtoolSetKernel(Mat,MatHtoolKernel,void*);
1908: PETSC_EXTERN PetscErrorCode MatHtoolGetPermutationSource(Mat,IS*);
1909: PETSC_EXTERN PetscErrorCode MatHtoolGetPermutationTarget(Mat,IS*);
1910: PETSC_EXTERN PetscErrorCode MatHtoolUsePermutation(Mat,PetscBool);
1911: /*E
1912:      MatHtoolCompressorType - Indicates the type of compressor used by a MATHTOOL

1914:    Level: beginner

1916:     Values:
1917: +   MAT_HTOOL_COMPRESSOR_SYMPARTIAL_ACA (default) - symmetric partial adaptive cross approximation
1918: .   MAT_HTOOL_COMPRESSOR_FULL_ACA - full adaptive cross approximation
1919: -   MAT_HTOOL_COMPRESSOR_SVD - singular value decomposition

1921: .seealso: MatCreateHtoolFromKernel(), MATHTOOL, MatHtoolClusteringType
1922: E*/
1923: typedef enum { MAT_HTOOL_COMPRESSOR_SYMPARTIAL_ACA, MAT_HTOOL_COMPRESSOR_FULL_ACA, MAT_HTOOL_COMPRESSOR_SVD } MatHtoolCompressorType;
1924: /*E
1925:      MatHtoolClusteringType - Indicates the type of clustering used by a MATHTOOL

1927:    Level: beginner

1929:     Values:
1930: +   MAT_HTOOL_CLUSTERING_PCA_REGULAR (default) - axis computed via principle component analysis, split uniformly
1931: .   MAT_HTOOL_CLUSTERING_PCA_GEOMETRIC - axis computed via principle component analysis, split barycentrically
1932: .   MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_REGULAR - axis along the largest extent of the bounding box, split uniformly
1933: -   MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_GEOMETRIC - axis along the largest extent of the bounding box, split barycentrically

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

1937: .seealso: MatCreateHtoolFromKernel(), MATHTOOL, MatHtoolCompressorType
1938: E*/
1939: 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;
1940: #endif

1942: /*
1943:    PETSc interface to MUMPS
1944: */
1945: #ifdef PETSC_HAVE_MUMPS
1946: PETSC_EXTERN PetscErrorCode MatMumpsSetIcntl(Mat,PetscInt,PetscInt);
1947: PETSC_EXTERN PetscErrorCode MatMumpsGetIcntl(Mat,PetscInt,PetscInt*);
1948: PETSC_EXTERN PetscErrorCode MatMumpsSetCntl(Mat,PetscInt,PetscReal);
1949: PETSC_EXTERN PetscErrorCode MatMumpsGetCntl(Mat,PetscInt,PetscReal*);

1951: PETSC_EXTERN PetscErrorCode MatMumpsGetInfo(Mat,PetscInt,PetscInt*);
1952: PETSC_EXTERN PetscErrorCode MatMumpsGetInfog(Mat,PetscInt,PetscInt*);
1953: PETSC_EXTERN PetscErrorCode MatMumpsGetRinfo(Mat,PetscInt,PetscReal*);
1954: PETSC_EXTERN PetscErrorCode MatMumpsGetRinfog(Mat,PetscInt,PetscReal*);
1955: PETSC_EXTERN PetscErrorCode MatMumpsGetInverse(Mat,Mat);
1956: PETSC_EXTERN PetscErrorCode MatMumpsGetInverseTranspose(Mat,Mat);
1957: #endif

1959: /*
1960:    PETSc interface to Mkl_Pardiso
1961: */
1962: #ifdef PETSC_HAVE_MKL_PARDISO
1963: PETSC_EXTERN PetscErrorCode MatMkl_PardisoSetCntl(Mat,PetscInt,PetscInt);
1964: #endif

1966: /*
1967:    PETSc interface to Mkl_CPardiso
1968: */
1969: #ifdef PETSC_HAVE_MKL_CPARDISO
1970: PETSC_EXTERN PetscErrorCode MatMkl_CPardisoSetCntl(Mat,PetscInt,PetscInt);
1971: #endif

1973: /*
1974:    PETSc interface to SUPERLU
1975: */
1976: #ifdef PETSC_HAVE_SUPERLU
1977: PETSC_EXTERN PetscErrorCode MatSuperluSetILUDropTol(Mat,PetscReal);
1978: #endif

1980: /*
1981:    PETSc interface to SUPERLU_DIST
1982: */
1983: #ifdef PETSC_HAVE_SUPERLU_DIST
1984: PETSC_EXTERN PetscErrorCode MatSuperluDistGetDiagU(Mat,PetscScalar*);
1985: #endif

1987: /*
1988:    PETSc interface to STRUMPACK
1989: */
1990: #ifdef PETSC_HAVE_STRUMPACK
1991: /*E
1992:     MatSTRUMPACKReordering - sparsity reducing ordering to be used in STRUMPACK

1994:     Level: intermediate
1995: E*/
1996: typedef enum {MAT_STRUMPACK_NATURAL=0,
1997:               MAT_STRUMPACK_METIS=1,
1998:               MAT_STRUMPACK_PARMETIS=2,
1999:               MAT_STRUMPACK_SCOTCH=3,
2000:               MAT_STRUMPACK_PTSCOTCH=4,
2001:               MAT_STRUMPACK_RCM=5} MatSTRUMPACKReordering;
2002: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetReordering(Mat,MatSTRUMPACKReordering);
2003: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetColPerm(Mat,PetscBool);
2004: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSRelTol(Mat,PetscReal);
2005: PETSC_DEPRECATED_FUNCTION("Use MatSTRUMPACKSetHSSRelTol() (since version 3.9)") static inline PetscErrorCode MatSTRUMPACKSetHSSRelCompTol(Mat mat,PetscReal rtol) {return MatSTRUMPACKSetHSSRelTol(mat,rtol);}
2006: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSAbsTol(Mat,PetscReal);
2007: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSMinSepSize(Mat,PetscInt);
2008: PETSC_DEPRECATED_FUNCTION("Use MatSTRUMPACKSetHSSMinSepSize() (since version 3.9)") static inline PetscErrorCode MatSTRUMPACKSetHSSMinSize(Mat mat,PetscInt hssminsize) {return MatSTRUMPACKSetHSSMinSepSize(mat,hssminsize);}
2009: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSMaxRank(Mat,PetscInt);
2010: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSLeafSize(Mat,PetscInt);
2011: #endif

2013: PETSC_EXTERN PetscErrorCode MatBindToCPU(Mat,PetscBool);
2014: PETSC_EXTERN PetscErrorCode MatBoundToCPU(Mat,PetscBool*);
2015: PETSC_DEPRECATED_FUNCTION("Use MatBindToCPU (since v3.13)") static inline PetscErrorCode MatPinToCPU(Mat A,PetscBool flg) {return MatBindToCPU(A,flg);}
2016: PETSC_EXTERN PetscErrorCode MatSetBindingPropagates(Mat,PetscBool);
2017: PETSC_EXTERN PetscErrorCode MatGetBindingPropagates(Mat,PetscBool *);

2019: typedef struct _n_SplitCSRMat *PetscSplitCSRDataStructure;
2020: PETSC_EXTERN PetscErrorCode MatCUSPARSEGetDeviceMatWrite(Mat,PetscSplitCSRDataStructure*);

2022: #ifdef PETSC_HAVE_KOKKOS_KERNELS
2023: PETSC_EXTERN PetscErrorCode MatKokkosGetDeviceMatWrite(Mat,PetscSplitCSRDataStructure*);
2024: PETSC_EXTERN PetscErrorCode MatSeqAIJKokkosSetDeviceMat(Mat,PetscSplitCSRDataStructure);
2025: PETSC_EXTERN PetscErrorCode MatSeqAIJKokkosGetDeviceMat(Mat,PetscSplitCSRDataStructure*);
2026: #endif

2028: #ifdef PETSC_HAVE_CUDA
2029: /*E
2030:     MatCUSPARSEStorageFormat - indicates the storage format for CUSPARSE (GPU)
2031:     matrices.

2033:     Not Collective

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

2039:     Level: intermediate

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

2043: .seealso: MatCUSPARSESetFormat(), MatCUSPARSEFormatOperation
2044: E*/

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

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

2051: /*E
2052:     MatCUSPARSEFormatOperation - indicates the operation of CUSPARSE (GPU)
2053:     matrices whose operation should use a particular storage format.

2055:     Not Collective

2057: +   MAT_CUSPARSE_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult
2058: .   MAT_CUSPARSE_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult
2059: .   MAT_CUSPARSE_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult
2060: -   MAT_CUSPARSE_ALL - sets the storage format for all CUSPARSE (GPU) matrices

2062:     Level: intermediate

2064: .seealso: MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat
2065: E*/
2066: typedef enum {MAT_CUSPARSE_MULT_DIAG, MAT_CUSPARSE_MULT_OFFDIAG, MAT_CUSPARSE_MULT, MAT_CUSPARSE_ALL} MatCUSPARSEFormatOperation;

2068: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
2069: PETSC_EXTERN PetscErrorCode MatCreateAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
2070: PETSC_EXTERN PetscErrorCode MatCUSPARSESetFormat(Mat,MatCUSPARSEFormatOperation,MatCUSPARSEStorageFormat);
2071: PETSC_EXTERN PetscErrorCode MatCUSPARSESetUseCPUSolve(Mat,PetscBool);
2072: typedef struct Mat_SeqAIJCUSPARSETriFactors* Mat_SeqAIJCUSPARSETriFactors_p;
2073: PETSC_EXTERN PetscErrorCode MatSeqAIJCUSPARSEGetIJ(Mat,PetscBool,const int**,const int**);
2074: PETSC_EXTERN PetscErrorCode MatSeqAIJCUSPARSERestoreIJ(Mat,PetscBool,const int**,const int**);
2075: PETSC_EXTERN PetscErrorCode MatSeqAIJCUSPARSEGetArrayRead(Mat,const PetscScalar**);
2076: PETSC_EXTERN PetscErrorCode MatSeqAIJCUSPARSERestoreArrayRead(Mat,const PetscScalar**);
2077: PETSC_EXTERN PetscErrorCode MatSeqAIJCUSPARSEGetArrayWrite(Mat,PetscScalar**);
2078: PETSC_EXTERN PetscErrorCode MatSeqAIJCUSPARSERestoreArrayWrite(Mat,PetscScalar**);
2079: PETSC_EXTERN PetscErrorCode MatSeqAIJCUSPARSEGetArray(Mat,PetscScalar**);
2080: PETSC_EXTERN PetscErrorCode MatSeqAIJCUSPARSERestoreArray(Mat,PetscScalar**);

2082: PETSC_EXTERN PetscErrorCode MatCreateDenseCUDA(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
2083: PETSC_EXTERN PetscErrorCode MatCreateSeqDenseCUDA(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
2084: PETSC_EXTERN PetscErrorCode MatMPIDenseCUDASetPreallocation(Mat,PetscScalar[]);
2085: PETSC_EXTERN PetscErrorCode MatSeqDenseCUDASetPreallocation(Mat,PetscScalar[]);
2086: PETSC_EXTERN PetscErrorCode MatDenseCUDAGetArrayWrite(Mat,PetscScalar**);
2087: PETSC_EXTERN PetscErrorCode MatDenseCUDAGetArrayRead(Mat,const PetscScalar**);
2088: PETSC_EXTERN PetscErrorCode MatDenseCUDAGetArray(Mat,PetscScalar**);
2089: PETSC_EXTERN PetscErrorCode MatDenseCUDARestoreArrayWrite(Mat,PetscScalar**);
2090: PETSC_EXTERN PetscErrorCode MatDenseCUDARestoreArrayRead(Mat,const PetscScalar**);
2091: PETSC_EXTERN PetscErrorCode MatDenseCUDARestoreArray(Mat,PetscScalar**);
2092: PETSC_EXTERN PetscErrorCode MatDenseCUDAPlaceArray(Mat,const PetscScalar*);
2093: PETSC_EXTERN PetscErrorCode MatDenseCUDAReplaceArray(Mat,const PetscScalar*);
2094: PETSC_EXTERN PetscErrorCode MatDenseCUDAResetArray(Mat);

2096: #endif

2098: #if defined(PETSC_HAVE_VIENNACL)
2099: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
2100: PETSC_EXTERN PetscErrorCode MatCreateAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
2101: #endif

2103: /*
2104:    PETSc interface to FFTW
2105: */
2106: #if defined(PETSC_HAVE_FFTW)
2107: PETSC_EXTERN PetscErrorCode VecScatterPetscToFFTW(Mat,Vec,Vec);
2108: PETSC_EXTERN PetscErrorCode VecScatterFFTWToPetsc(Mat,Vec,Vec);
2109: PETSC_EXTERN PetscErrorCode MatCreateVecsFFTW(Mat,Vec*,Vec*,Vec*);
2110: #endif

2112: #if defined(PETSC_HAVE_SCALAPACK)
2113: PETSC_EXTERN PetscErrorCode MatCreateScaLAPACK(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
2114: PETSC_EXTERN PetscErrorCode MatScaLAPACKSetBlockSizes(Mat,PetscInt,PetscInt);
2115: PETSC_EXTERN PetscErrorCode MatScaLAPACKGetBlockSizes(Mat,PetscInt*,PetscInt*);
2116: #endif

2118: PETSC_EXTERN PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*);
2119: PETSC_EXTERN PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*);
2120: PETSC_EXTERN PetscErrorCode MatNestGetISs(Mat,IS[],IS[]);
2121: PETSC_EXTERN PetscErrorCode MatNestGetLocalISs(Mat,IS[],IS[]);
2122: PETSC_EXTERN PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***);
2123: PETSC_EXTERN PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*);
2124: PETSC_EXTERN PetscErrorCode MatNestSetVecType(Mat,VecType);
2125: PETSC_EXTERN PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]);
2126: PETSC_EXTERN PetscErrorCode MatNestSetSubMat(Mat,PetscInt,PetscInt,Mat);

2128: PETSC_EXTERN PetscErrorCode MatChop(Mat,PetscReal);
2129: PETSC_EXTERN PetscErrorCode MatComputeBandwidth(Mat,PetscReal,PetscInt*);

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

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

2135: PETSC_INTERN PetscErrorCode MatHeaderMerge(Mat,Mat*);
2136: PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat,Mat*);

2138: PETSC_EXTERN PetscErrorCode MatSeqAIJGetCSRAndMemType(Mat,const PetscInt**,const PetscInt**,PetscScalar**,PetscMemType*);
2139: #endif