Actual source code: petscmat.h

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

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

 12:    Level: beginner

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

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

 21:    Level: beginner

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

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

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

118:    Level: beginner

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

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

146: /*E
147:     MatFactorType - indicates what type of factorization is requested

149:     Level: beginner

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

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

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

173: /*E
174:     MatProductType - indicates what type of matrix product is requested

176:     Level: beginner

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

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

186:    Level: beginner

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

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

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

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

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

225:     Level: beginner

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

229: .seealso: MatCreateSubMatrices(), MatCreateSubMatrix(), MatDestroyMatrices(), MatConvert()
230: E*/
231: typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX,MAT_IGNORE_MATRIX,MAT_INPLACE_MATRIX} MatReuse;

233: /*E
234:     MatCreateSubMatrixOption - Indicates if matrices obtained from a call to MatCreateSubMatrices()
235:      include the matrix values. Currently it is only used by MatGetSeqNonzerostructure().

237:     Level: beginner

239: .seealso: MatGetSeqNonzerostructure()
240: E*/
241: typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatCreateSubMatrixOption;

243: PETSC_EXTERN PetscErrorCode MatInitializePackage(void);

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

259: PETSC_EXTERN PetscFunctionList MatList;
260: PETSC_EXTERN PetscFunctionList MatColoringList;
261: PETSC_EXTERN PetscFunctionList MatPartitioningList;

263: /*E
264:     MatStructure - Indicates if two matrices have the same nonzero structure

266:     Level: beginner

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

270: .seealso: MatCopy(), MatAXPY()
271: E*/
272: typedef enum {DIFFERENT_NONZERO_PATTERN,SUBSET_NONZERO_PATTERN,SAME_NONZERO_PATTERN} MatStructure;

274: #if defined PETSC_HAVE_MKL_SPARSE
275: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJMKL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
276: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJMKL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
277: PETSC_EXTERN PetscErrorCode MatCreateBAIJMKL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
278: PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJMKL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
279: #endif

281: PETSC_EXTERN PetscErrorCode MatCreateSeqSELL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
282: PETSC_EXTERN PetscErrorCode MatCreateSELL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
283: PETSC_EXTERN PetscErrorCode MatSeqSELLSetPreallocation(Mat,PetscInt,const PetscInt[]);
284: PETSC_EXTERN PetscErrorCode MatMPISELLSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);

286: PETSC_EXTERN PetscErrorCode MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
287: PETSC_EXTERN PetscErrorCode MatCreateDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
288: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
289: PETSC_EXTERN PetscErrorCode MatCreateAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
290: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
291: PETSC_EXTERN PetscErrorCode MatUpdateMPIAIJWithArrays(Mat,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
292: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*);
293: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithSeqAIJ(MPI_Comm,Mat,Mat,const PetscInt[],Mat*);

295: PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
296: PETSC_EXTERN PetscErrorCode MatCreateBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
297: PETSC_EXTERN PetscErrorCode MatCreateMPIBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat*);

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

302: PETSC_EXTERN PetscErrorCode MatCreateSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
303: PETSC_EXTERN PetscErrorCode MatCreateMPISBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
304: PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
305: PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
306: PETSC_EXTERN PetscErrorCode MatXAIJSetPreallocation(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscInt[],const PetscInt[]);

308: PETSC_EXTERN PetscErrorCode MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*);
309: PETSC_EXTERN PetscErrorCode MatCreateNormal(Mat,Mat*);
310: PETSC_EXTERN PetscErrorCode MatCreateNormalHermitian(Mat,Mat*);
311: PETSC_EXTERN PetscErrorCode MatCreateLRC(Mat,Mat,Vec,Mat,Mat*);
312: PETSC_EXTERN PetscErrorCode MatLRCGetMats(Mat,Mat*,Mat*,Vec*,Mat*);
313: PETSC_EXTERN PetscErrorCode MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,ISLocalToGlobalMapping,Mat*);
314: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
315: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);

317: PETSC_EXTERN PetscErrorCode MatCreateScatter(MPI_Comm,VecScatter,Mat*);
318: PETSC_EXTERN PetscErrorCode MatScatterSetVecScatter(Mat,VecScatter);
319: PETSC_EXTERN PetscErrorCode MatScatterGetVecScatter(Mat,VecScatter*);
320: PETSC_EXTERN PetscErrorCode MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*);
321: PETSC_EXTERN PetscErrorCode MatCompositeAddMat(Mat,Mat);
322: PETSC_EXTERN PetscErrorCode MatCompositeMerge(Mat);
323: typedef enum {MAT_COMPOSITE_MERGE_RIGHT,MAT_COMPOSITE_MERGE_LEFT} MatCompositeMergeType;
324: PETSC_EXTERN PetscErrorCode MatCompositeSetMergeType(Mat,MatCompositeMergeType);
325: PETSC_EXTERN PetscErrorCode MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*);
326: typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType;
327: PETSC_EXTERN PetscErrorCode MatCompositeSetType(Mat,MatCompositeType);
328: PETSC_EXTERN PetscErrorCode MatCompositeGetType(Mat,MatCompositeType*);
329: PETSC_EXTERN PetscErrorCode MatCompositeSetMatStructure(Mat,MatStructure);
330: PETSC_EXTERN PetscErrorCode MatCompositeGetMatStructure(Mat,MatStructure*);
331: PETSC_EXTERN PetscErrorCode MatCompositeGetNumberMat(Mat,PetscInt*);
332: PETSC_EXTERN PetscErrorCode MatCompositeGetMat(Mat,PetscInt,Mat*);
333: PETSC_EXTERN PetscErrorCode MatCompositeSetScalings(Mat,const PetscScalar*);

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

338: PETSC_EXTERN PetscErrorCode MatCreateTranspose(Mat,Mat*);
339: PETSC_EXTERN PetscErrorCode MatTransposeGetMat(Mat,Mat*);
340: PETSC_EXTERN PetscErrorCode MatCreateHermitianTranspose(Mat,Mat*);
341: PETSC_EXTERN PetscErrorCode MatHermitianTransposeGetMat(Mat,Mat*);
342: PETSC_EXTERN PetscErrorCode MatCreateSubMatrixVirtual(Mat,IS,IS,Mat*);
343: PETSC_EXTERN PetscErrorCode MatSubMatrixVirtualUpdate(Mat,Mat,IS,IS);
344: PETSC_EXTERN PetscErrorCode MatCreateLocalRef(Mat,IS,IS,Mat*);
345: PETSC_EXTERN PetscErrorCode MatCreateConstantDiagonal(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar,Mat*);

347: #if defined(PETSC_HAVE_HYPRE)
348: PETSC_EXTERN PetscErrorCode MatHYPRESetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
349: #endif

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

353: PETSC_EXTERN PetscErrorCode MatResetPreallocation(Mat);
354: PETSC_EXTERN PetscErrorCode MatSetUp(Mat);
355: PETSC_EXTERN PetscErrorCode MatDestroy(Mat*);
356: PETSC_EXTERN PetscErrorCode MatGetNonzeroState(Mat,PetscObjectState*);

358: PETSC_EXTERN PetscErrorCode MatConjugate(Mat);
359: PETSC_EXTERN PetscErrorCode MatRealPart(Mat);
360: PETSC_EXTERN PetscErrorCode MatImaginaryPart(Mat);
361: PETSC_EXTERN PetscErrorCode MatGetDiagonalBlock(Mat,Mat*);
362: PETSC_EXTERN PetscErrorCode MatGetTrace(Mat,PetscScalar*);
363: PETSC_EXTERN PetscErrorCode MatInvertBlockDiagonal(Mat,const PetscScalar **);
364: PETSC_EXTERN PetscErrorCode MatInvertVariableBlockDiagonal(Mat,PetscInt,const PetscInt*,PetscScalar*);
365: PETSC_EXTERN PetscErrorCode MatInvertBlockDiagonalMat(Mat,Mat);

367: /* ------------------------------------------------------------*/
368: PETSC_EXTERN PetscErrorCode MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
369: PETSC_EXTERN PetscErrorCode MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
370: PETSC_EXTERN PetscErrorCode MatSetValuesRow(Mat,PetscInt,const PetscScalar[]);
371: PETSC_EXTERN PetscErrorCode MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]);
372: PETSC_EXTERN PetscErrorCode MatSetValuesBatch(Mat,PetscInt,PetscInt,PetscInt[],const PetscScalar[]);
373: PETSC_EXTERN PetscErrorCode MatSetRandom(Mat,PetscRandom);

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

379:    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).
380:    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.

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

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

386:    Level: beginner

388: .seealso:  MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockedStencil(), DMDAVecGetArray(), DMDAVecGetArrayF90()
389: S*/
390: typedef struct {
391:   PetscInt k,j,i,c;
392: } MatStencil;

394: PETSC_EXTERN PetscErrorCode MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
395: PETSC_EXTERN PetscErrorCode MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
396: PETSC_EXTERN PetscErrorCode MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt);

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

402:     Level: beginner

404: .seealso: MatAssemblyBegin(), MatAssemblyEnd()
405: E*/
406: typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType;
407: PETSC_EXTERN PetscErrorCode MatAssemblyBegin(Mat,MatAssemblyType);
408: PETSC_EXTERN PetscErrorCode MatAssemblyEnd(Mat,MatAssemblyType);
409: PETSC_EXTERN PetscErrorCode MatAssembled(Mat,PetscBool *);



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

416:     Level: beginner

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

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

424: .seealso: MatSetOption()
425: E*/
426: typedef enum {MAT_OPTION_MIN = -3,
427:               MAT_UNUSED_NONZERO_LOCATION_ERR = -2,
428:               MAT_ROW_ORIENTED = -1,
429:               MAT_SYMMETRIC = 1,
430:               MAT_STRUCTURALLY_SYMMETRIC = 2,
431:               MAT_NEW_DIAGONALS = 3,
432:               MAT_IGNORE_OFF_PROC_ENTRIES = 4,
433:               MAT_USE_HASH_TABLE = 5,
434:               MAT_KEEP_NONZERO_PATTERN = 6,
435:               MAT_IGNORE_ZERO_ENTRIES = 7,
436:               MAT_USE_INODES = 8,
437:               MAT_HERMITIAN = 9,
438:               MAT_SYMMETRY_ETERNAL = 10,
439:               MAT_NEW_NONZERO_LOCATION_ERR = 11,
440:               MAT_IGNORE_LOWER_TRIANGULAR = 12,
441:               MAT_ERROR_LOWER_TRIANGULAR = 13,
442:               MAT_GETROW_UPPERTRIANGULAR = 14,
443:               MAT_SPD = 15,
444:               MAT_NO_OFF_PROC_ZERO_ROWS = 16,
445:               MAT_NO_OFF_PROC_ENTRIES = 17,
446:               MAT_NEW_NONZERO_LOCATIONS = 18,
447:               MAT_NEW_NONZERO_ALLOCATION_ERR = 19,
448:               MAT_SUBSET_OFF_PROC_ENTRIES = 20,
449:               MAT_SUBMAT_SINGLEIS = 21,
450:               MAT_STRUCTURE_ONLY = 22,
451:               MAT_SORTED_FULL = 23,
452:               MAT_OPTION_MAX = 24} MatOption;

454: PETSC_EXTERN const char *const *MatOptions;
455: PETSC_EXTERN PetscErrorCode MatSetOption(Mat,MatOption,PetscBool);
456: PETSC_EXTERN PetscErrorCode MatGetOption(Mat,MatOption,PetscBool*);
457: PETSC_EXTERN PetscErrorCode MatPropagateSymmetryOptions(Mat,Mat);
458: PETSC_EXTERN PetscErrorCode MatGetType(Mat,MatType*);

460: PETSC_EXTERN PetscErrorCode MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
461: PETSC_EXTERN PetscErrorCode MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
462: PETSC_EXTERN PetscErrorCode MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
463: PETSC_EXTERN PetscErrorCode MatGetRowUpperTriangular(Mat);
464: PETSC_EXTERN PetscErrorCode MatRestoreRowUpperTriangular(Mat);
465: PETSC_EXTERN PetscErrorCode MatGetColumnVector(Mat,Vec,PetscInt);
466: PETSC_EXTERN PetscErrorCode MatSeqAIJGetArray(Mat,PetscScalar *[]);
467: PETSC_EXTERN PetscErrorCode MatSeqAIJGetArrayRead(Mat,const PetscScalar *[]);
468: PETSC_EXTERN PetscErrorCode MatSeqAIJRestoreArray(Mat,PetscScalar *[]);
469: PETSC_EXTERN PetscErrorCode MatSeqAIJRestoreArrayRead(Mat,const PetscScalar *[]);
470: PETSC_EXTERN PetscErrorCode MatSeqAIJGetMaxRowNonzeros(Mat,PetscInt*);
471: PETSC_EXTERN PetscErrorCode MatSeqAIJSetValuesLocalFast(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
472: PETSC_EXTERN PetscErrorCode MatSeqAIJSetType(Mat,MatType);
473: PETSC_EXTERN PetscErrorCode MatSeqAIJRegister(const char[],PetscErrorCode (*)(Mat,MatType,MatReuse,Mat *));
474: PETSC_EXTERN PetscFunctionList MatSeqAIJList;
475: PETSC_EXTERN PetscErrorCode MatSeqBAIJGetArray(Mat,PetscScalar *[]);
476: PETSC_EXTERN PetscErrorCode MatSeqBAIJRestoreArray(Mat,PetscScalar *[]);
477: PETSC_EXTERN PetscErrorCode MatSeqSBAIJGetArray(Mat,PetscScalar *[]);
478: PETSC_EXTERN PetscErrorCode MatSeqSBAIJRestoreArray(Mat,PetscScalar *[]);
479: PETSC_EXTERN PetscErrorCode MatDenseGetArray(Mat,PetscScalar *[]);
480: PETSC_EXTERN PetscErrorCode MatDenseRestoreArray(Mat,PetscScalar *[]);
481: PETSC_EXTERN PetscErrorCode MatDensePlaceArray(Mat,const PetscScalar[]);
482: PETSC_EXTERN PetscErrorCode MatDenseReplaceArray(Mat,const PetscScalar[]);
483: PETSC_EXTERN PetscErrorCode MatDenseResetArray(Mat);
484: PETSC_EXTERN PetscErrorCode MatDenseGetArrayRead(Mat,const PetscScalar *[]);
485: PETSC_EXTERN PetscErrorCode MatDenseRestoreArrayRead(Mat,const PetscScalar *[]);
486: PETSC_EXTERN PetscErrorCode MatDenseGetArrayWrite(Mat,PetscScalar *[]);
487: PETSC_EXTERN PetscErrorCode MatDenseRestoreArrayWrite(Mat,PetscScalar *[]);
488: PETSC_EXTERN PetscErrorCode MatGetBlockSize(Mat,PetscInt *);
489: PETSC_EXTERN PetscErrorCode MatSetBlockSize(Mat,PetscInt);
490: PETSC_EXTERN PetscErrorCode MatGetBlockSizes(Mat,PetscInt *,PetscInt *);
491: PETSC_EXTERN PetscErrorCode MatSetBlockSizes(Mat,PetscInt,PetscInt);
492: PETSC_EXTERN PetscErrorCode MatSetBlockSizesFromMats(Mat,Mat,Mat);
493: PETSC_EXTERN PetscErrorCode MatSetVariableBlockSizes(Mat,PetscInt,PetscInt*);
494: PETSC_EXTERN PetscErrorCode MatGetVariableBlockSizes(Mat,PetscInt*,const PetscInt**);

496: PETSC_EXTERN PetscErrorCode MatDenseGetColumn(Mat,PetscInt,PetscScalar*[]);
497: PETSC_EXTERN PetscErrorCode MatDenseRestoreColumn(Mat,PetscScalar*[]);
498: PETSC_EXTERN PetscErrorCode MatDenseGetColumnVec(Mat,PetscInt,Vec*);
499: PETSC_EXTERN PetscErrorCode MatDenseRestoreColumnVec(Mat,PetscInt,Vec*);
500: PETSC_EXTERN PetscErrorCode MatDenseGetColumnVecRead(Mat,PetscInt,Vec*);
501: PETSC_EXTERN PetscErrorCode MatDenseRestoreColumnVecRead(Mat,PetscInt,Vec*);
502: PETSC_EXTERN PetscErrorCode MatDenseGetColumnVecWrite(Mat,PetscInt,Vec*);
503: PETSC_EXTERN PetscErrorCode MatDenseRestoreColumnVecWrite(Mat,PetscInt,Vec*);
504: PETSC_EXTERN PetscErrorCode MatDenseGetSubMatrix(Mat,PetscInt,PetscInt,Mat*);
505: PETSC_EXTERN PetscErrorCode MatDenseRestoreSubMatrix(Mat,Mat*);

507: PETSC_EXTERN PetscErrorCode MatMult(Mat,Vec,Vec);
508: PETSC_EXTERN PetscErrorCode MatMultDiagonalBlock(Mat,Vec,Vec);
509: PETSC_EXTERN PetscErrorCode MatMultAdd(Mat,Vec,Vec,Vec);
510: PETSC_EXTERN PetscErrorCode MatMultTranspose(Mat,Vec,Vec);
511: PETSC_EXTERN PetscErrorCode MatMultHermitianTranspose(Mat,Vec,Vec);
512: PETSC_EXTERN PetscErrorCode MatIsTranspose(Mat,Mat,PetscReal,PetscBool *);
513: PETSC_EXTERN PetscErrorCode MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *);
514: PETSC_EXTERN PetscErrorCode MatMultTransposeAdd(Mat,Vec,Vec,Vec);
515: PETSC_EXTERN PetscErrorCode MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec);
516: PETSC_EXTERN PetscErrorCode MatMultConstrained(Mat,Vec,Vec);
517: PETSC_EXTERN PetscErrorCode MatMultTransposeConstrained(Mat,Vec,Vec);
518: PETSC_EXTERN PetscErrorCode MatMatSolve(Mat,Mat,Mat);
519: PETSC_EXTERN PetscErrorCode MatMatSolveTranspose(Mat,Mat,Mat);
520: PETSC_EXTERN PetscErrorCode MatMatTransposeSolve(Mat,Mat,Mat);
521: PETSC_EXTERN PetscErrorCode MatResidual(Mat,Vec,Vec,Vec);

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

527:     Level: beginner

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

531: $   MAT_DO_NOT_COPY_VALUES    - Create a matrix using the same nonzero pattern as the original matrix,
532: $                               with zeros for the numerical values.
533: $   MAT_COPY_VALUES           - Create a matrix with the same nonzero pattern as the original matrix
534: $                               and with the same numerical values.
535: $   MAT_SHARE_NONZERO_PATTERN - Create a matrix that shares the nonzero structure with the previous matrix
536: $                               and does not copy it, using zeros for the numerical values. The parent and
537: $                               child matrices will share their index (i and j) arrays, and you cannot
538: $                               insert new nonzero entries into either matrix.

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

544: .seealso: MatDuplicate()
545: E*/
546: typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption;

548: PETSC_EXTERN PetscErrorCode MatConvert(Mat,MatType,MatReuse,Mat*);
549: PETSC_EXTERN PetscErrorCode MatDuplicate(Mat,MatDuplicateOption,Mat*);


552: PETSC_EXTERN PetscErrorCode MatCopy(Mat,Mat,MatStructure);
553: PETSC_EXTERN PetscErrorCode MatView(Mat,PetscViewer);
554: PETSC_EXTERN PetscErrorCode MatIsSymmetric(Mat,PetscReal,PetscBool *);
555: PETSC_EXTERN PetscErrorCode MatIsStructurallySymmetric(Mat,PetscBool *);
556: PETSC_EXTERN PetscErrorCode MatIsHermitian(Mat,PetscReal,PetscBool *);
557: PETSC_EXTERN PetscErrorCode MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *);
558: PETSC_EXTERN PetscErrorCode MatIsHermitianKnown(Mat,PetscBool *,PetscBool *);
559: PETSC_EXTERN PetscErrorCode MatMissingDiagonal(Mat,PetscBool  *,PetscInt *);
560: PETSC_EXTERN PetscErrorCode MatLoad(Mat, PetscViewer);

562: PETSC_EXTERN PetscErrorCode MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
563: PETSC_EXTERN PetscErrorCode MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool  *);
564: PETSC_EXTERN PetscErrorCode MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
565: PETSC_EXTERN PetscErrorCode MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool  *);

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

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

572:    Level: intermediate

574: .seealso:  MatGetInfo(), MatInfoType
575: S*/
576: typedef struct {
577:   PetscLogDouble block_size;                         /* block size */
578:   PetscLogDouble nz_allocated,nz_used,nz_unneeded;   /* number of nonzeros */
579:   PetscLogDouble memory;                             /* memory allocated */
580:   PetscLogDouble assemblies;                         /* number of matrix assemblies called */
581:   PetscLogDouble mallocs;                            /* number of mallocs during MatSetValues() */
582:   PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
583:   PetscLogDouble factor_mallocs;                     /* number of mallocs during factorization */
584: } MatInfo;

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

590:     Level: beginner

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

594: .seealso: MatGetInfo(), MatInfo
595: E*/
596: typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
597: PETSC_EXTERN PetscErrorCode MatGetInfo(Mat,MatInfoType,MatInfo*);
598: PETSC_EXTERN PetscErrorCode MatGetDiagonal(Mat,Vec);
599: PETSC_EXTERN PetscErrorCode MatGetRowMax(Mat,Vec,PetscInt[]);
600: PETSC_EXTERN PetscErrorCode MatGetRowMin(Mat,Vec,PetscInt[]);
601: PETSC_EXTERN PetscErrorCode MatGetRowMaxAbs(Mat,Vec,PetscInt[]);
602: PETSC_EXTERN PetscErrorCode MatGetRowMinAbs(Mat,Vec,PetscInt[]);
603: PETSC_EXTERN PetscErrorCode MatGetRowSum(Mat,Vec);
604: PETSC_EXTERN PetscErrorCode MatTranspose(Mat,MatReuse,Mat*);
605: PETSC_EXTERN PetscErrorCode MatHermitianTranspose(Mat,MatReuse,Mat*);
606: PETSC_EXTERN PetscErrorCode MatPermute(Mat,IS,IS,Mat*);
607: PETSC_EXTERN PetscErrorCode MatDiagonalScale(Mat,Vec,Vec);
608: PETSC_EXTERN PetscErrorCode MatDiagonalSet(Mat,Vec,InsertMode);

610: PETSC_EXTERN PetscErrorCode MatEqual(Mat,Mat,PetscBool*);
611: PETSC_EXTERN PetscErrorCode MatMultEqual(Mat,Mat,PetscInt,PetscBool*);
612: PETSC_EXTERN PetscErrorCode MatMultAddEqual(Mat,Mat,PetscInt,PetscBool*);
613: PETSC_EXTERN PetscErrorCode MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool*);
614: PETSC_EXTERN PetscErrorCode MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool*);
615: PETSC_EXTERN PetscErrorCode MatMatMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*);
616: PETSC_EXTERN PetscErrorCode MatTransposeMatMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*);
617: PETSC_EXTERN PetscErrorCode MatMatTransposeMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*);
618: PETSC_EXTERN PetscErrorCode MatPtAPMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*);
619: PETSC_EXTERN PetscErrorCode MatRARtMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*);
620: PETSC_EXTERN PetscErrorCode MatIsLinear(Mat,PetscInt,PetscBool*);

622: PETSC_EXTERN PetscErrorCode MatNorm(Mat,NormType,PetscReal*);
623: PETSC_EXTERN PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal*);
624: PETSC_EXTERN PetscErrorCode MatZeroEntries(Mat);
625: PETSC_EXTERN PetscErrorCode MatSetInf(Mat);
626: PETSC_EXTERN PetscErrorCode MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
627: PETSC_EXTERN PetscErrorCode MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec);
628: PETSC_EXTERN PetscErrorCode MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec);
629: PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsStencil(Mat,PetscInt,const MatStencil[],PetscScalar,Vec,Vec);
630: PETSC_EXTERN PetscErrorCode MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
631: PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec);

633: PETSC_EXTERN PetscErrorCode MatGetSize(Mat,PetscInt*,PetscInt*);
634: PETSC_EXTERN PetscErrorCode MatGetLocalSize(Mat,PetscInt*,PetscInt*);
635: PETSC_EXTERN PetscErrorCode MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
636: PETSC_EXTERN PetscErrorCode MatGetOwnershipRanges(Mat,const PetscInt**);
637: PETSC_EXTERN PetscErrorCode MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*);
638: PETSC_EXTERN PetscErrorCode MatGetOwnershipRangesColumn(Mat,const PetscInt**);
639: PETSC_EXTERN PetscErrorCode MatGetOwnershipIS(Mat,IS*,IS*);

641: PETSC_EXTERN PetscErrorCode MatCreateSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
642: PETSC_DEPRECATED_FUNCTION("Use MatCreateSubMatrices() (since version 3.8)") PETSC_STATIC_INLINE PetscErrorCode MatGetSubMatrices(Mat mat,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *submat[]) {return MatCreateSubMatrices(mat,n,irow,icol,scall,submat);}
643: PETSC_EXTERN PetscErrorCode MatCreateSubMatricesMPI(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
644: PETSC_DEPRECATED_FUNCTION("Use MatCreateSubMatricesMPI() (since version 3.8)") PETSC_STATIC_INLINE PetscErrorCode MatGetSubMatricesMPI(Mat mat,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *submat[]) {return MatCreateSubMatricesMPI(mat,n,irow,icol,scall,submat);}
645: PETSC_EXTERN PetscErrorCode MatDestroyMatrices(PetscInt,Mat *[]);
646: PETSC_EXTERN PetscErrorCode MatDestroySubMatrices(PetscInt,Mat *[]);
647: PETSC_EXTERN PetscErrorCode MatCreateSubMatrix(Mat,IS,IS,MatReuse,Mat *);
648: PETSC_DEPRECATED_FUNCTION("Use MatCreateSubMatrix() (since version 3.8)") PETSC_STATIC_INLINE PetscErrorCode MatGetSubMatrix(Mat mat,IS isrow,IS iscol,MatReuse cll,Mat *newmat) {return MatCreateSubMatrix(mat,isrow,iscol,cll,newmat);}
649: PETSC_EXTERN PetscErrorCode MatGetLocalSubMatrix(Mat,IS,IS,Mat*);
650: PETSC_EXTERN PetscErrorCode MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*);
651: PETSC_EXTERN PetscErrorCode MatGetSeqNonzeroStructure(Mat,Mat*);
652: PETSC_EXTERN PetscErrorCode MatDestroySeqNonzeroStructure(Mat*);

654: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJ(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*);
655: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*);
656: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJNumeric(Mat,Mat);
657: PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMat(Mat,MatReuse,Mat*);
658: PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*);
659: PETSC_EXTERN PetscErrorCode MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,Mat*);
660: PETSC_EXTERN PetscErrorCode MatGetGhosts(Mat, PetscInt *,const PetscInt *[]);

662: PETSC_EXTERN PetscErrorCode MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt);
663: PETSC_EXTERN PetscErrorCode MatIncreaseOverlapSplit(Mat mat,PetscInt n,IS is[],PetscInt ov);
664: PETSC_EXTERN PetscErrorCode MatMPIAIJSetUseScalableIncreaseOverlap(Mat,PetscBool);

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

668: PETSC_EXTERN PetscErrorCode MatMatMatMult(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
669: PETSC_EXTERN PetscErrorCode MatGalerkin(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);

671: PETSC_EXTERN PetscErrorCode MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*);
672: PETSC_EXTERN PetscErrorCode MatRARt(Mat,Mat,MatReuse,PetscReal,Mat*);

674: PETSC_EXTERN PetscErrorCode MatTransposeMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
675: PETSC_EXTERN PetscErrorCode MatMatTransposeMult(Mat,Mat,MatReuse,PetscReal,Mat*);

677: PETSC_EXTERN PetscErrorCode MatAXPY(Mat,PetscScalar,Mat,MatStructure);
678: PETSC_EXTERN PetscErrorCode MatAYPX(Mat,PetscScalar,Mat,MatStructure);

680: PETSC_EXTERN PetscErrorCode MatScale(Mat,PetscScalar);
681: PETSC_EXTERN PetscErrorCode MatShift(Mat,PetscScalar);

683: PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
684: PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*);
685: PETSC_EXTERN PetscErrorCode MatGetLayouts(Mat,PetscLayout*,PetscLayout*);
686: PETSC_EXTERN PetscErrorCode MatSetLayouts(Mat,PetscLayout,PetscLayout);
687: PETSC_EXTERN PetscErrorCode MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
688: PETSC_EXTERN PetscErrorCode MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
689: PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
690: PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
691: PETSC_EXTERN PetscErrorCode MatGetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
692: PETSC_EXTERN PetscErrorCode MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
693: PETSC_EXTERN PetscErrorCode MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);

695: PETSC_EXTERN PetscErrorCode MatStashSetInitialSize(Mat,PetscInt,PetscInt);
696: PETSC_EXTERN PetscErrorCode MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);

698: PETSC_EXTERN PetscErrorCode MatInterpolate(Mat,Vec,Vec);
699: PETSC_EXTERN PetscErrorCode MatInterpolateAdd(Mat,Vec,Vec,Vec);
700: PETSC_EXTERN PetscErrorCode MatRestrict(Mat,Vec,Vec);
701: PETSC_EXTERN PetscErrorCode MatCreateVecs(Mat,Vec*,Vec*);
702: PETSC_DEPRECATED_FUNCTION("Use MatCreateVecs() (since version 3.6)") PETSC_STATIC_INLINE PetscErrorCode MatGetVecs(Mat mat,Vec *x,Vec *y) {return MatCreateVecs(mat,x,y);}
703: PETSC_EXTERN PetscErrorCode MatCreateRedundantMatrix(Mat,PetscInt,MPI_Comm,MatReuse,Mat*);
704: PETSC_EXTERN PetscErrorCode MatGetMultiProcBlock(Mat,MPI_Comm,MatReuse,Mat*);
705: PETSC_EXTERN PetscErrorCode MatFindZeroDiagonals(Mat,IS*);
706: PETSC_EXTERN PetscErrorCode MatFindOffBlockDiagonalEntries(Mat,IS*);
707: PETSC_EXTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);

709: /*MC
710:    MatSetValue - Set a single entry into a matrix.

712:    Not collective

714:    Synopsis:
715: #include <petscmat.h>
716:      PetscErrorCode MatSetValue(Mat m,PetscInt row,PetscInt col,PetscScalar value,InsertMode mode)

718:    Input Parameters:
719: +  m - the matrix
720: .  row - the row location of the entry
721: .  col - the column location of the entry
722: .  value - the value to insert
723: -  mode - either INSERT_VALUES or ADD_VALUES

725:    Notes:
726:    For efficiency one should use MatSetValues() and set several or many
727:    values simultaneously if possible.

729:    Level: beginner

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

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

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

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

743:    Synopsis:
744: #include <petscmat.h>
745:    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)

747:    Collective

749:    Input Parameters:
750: +  comm - the communicator that will share the eventually allocated matrix
751: .  nrows - the number of LOCAL rows in the matrix
752: -  ncols - the number of LOCAL columns in the matrix

754:    Output Parameters:
755: +  dnz - the array that will be passed to the matrix preallocation routines
756: -  onz - the other array passed to the matrix preallocation routines

758:    Level: intermediate

760:    Notes:
761:     See Users-Manual: ch_performance for more details.

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

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

767: .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
768:           MatPreallocateSymmetricSetLocalBlock()
769: M*/
770: #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
771: do { \
772:   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ncols = (ncols),__rstart,__start,__end; \
773:   _4_PetscCalloc2((size_t)__nrows,&dnz,(size_t)__nrows,&onz);CHKERRQ(_4_ierr); \
774:   __start = 0; __end = __start;                                         \
775:   _4_MPI_Scan(&__ncols,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ncols;\
776:   _4_MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows; \
777:   do { } while (0)

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

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

787:    Not Collective

789:    Input Parameters:
790: +  map - the row mapping from local numbering to global numbering
791: .  nrows - the number of rows indicated
792: .  rows - the indices of the rows
793: .  cmap - the column mapping from local to global numbering
794: .  ncols - the number of columns in the matrix
795: .  cols - the columns indicated
796: .  dnz - the array that will be passed to the matrix preallocation routines
797: -  onz - the other array passed to the matrix preallocation routines

799:    Level: intermediate

801:    Notes:
802:     See Users-Manual: ch_performance for more details.

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

806: .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock()
807:           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocalRemoveDups()
808: M*/
809: #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
810: do {\
811:   PetscInt __l;\
812:   _4_ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
813:   _4_ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
814:   for (__l=0;__l<nrows;__l++) {\
815:     _4_MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
816:   }\
817:  } while (0)

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

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

827:    Not Collective

829:    Input Parameters:
830: +  map - the row mapping from local numbering to global numbering
831: .  nrows - the number of rows indicated
832: .  rows - the indices of the rows (these values are mapped to the global values)
833: .  cmap - the column mapping from local to global numbering
834: .  ncols - the number of columns in the matrix   (this value will be changed if duplicate columns are found)
835: .  cols - the columns indicated (these values are mapped to the global values, they are then sorted and duplicates removed)
836: .  dnz - the array that will be passed to the matrix preallocation routines
837: -  onz - the other array passed to the matrix preallocation routines

839:    Level: intermediate

841:    Notes:
842:     See Users-Manual: ch_performance for more details.

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

846: .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock()
847:           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal()
848: M*/
849: #define MatPreallocateSetLocalRemoveDups(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
850: do {\
851:   PetscInt __l;\
852:   _4_ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
853:   _4_ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
854:   _4_PetscSortRemoveDupsInt(&ncols,cols);CHKERRQ(_4_ierr);\
855:   for (__l=0;__l<nrows;__l++) {\
856:     _4_MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
857:   }\
858:  } while (0)

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

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

868:    Not Collective

870:    Input Parameters:
871: +  map - the row mapping from local numbering to global numbering
872: .  nrows - the number of rows indicated
873: .  rows - the indices of the rows
874: .  cmap - the column mapping from local to global numbering
875: .  ncols - the number of columns in the matrix
876: .  cols - the columns indicated
877: .  dnz - the array that will be passed to the matrix preallocation routines
878: -  onz - the other array passed to the matrix preallocation routines

880:    Level: intermediate

882:    Notes:
883:     See Users-Manual: ch_performance for more details.

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

887: .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock()
888:           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock()
889: M*/
890: #define MatPreallocateSetLocalBlock(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
891: do {\
892:   PetscInt __l;\
893:   _4_ISLocalToGlobalMappingApplyBlock(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
894:   _4_ISLocalToGlobalMappingApplyBlock(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
895:   for (__l=0;__l<nrows;__l++) {\
896:     _4_MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
897:   }\
898: } while (0)

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

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

908:    Not Collective

910:    Input Parameters:
911: +  map - the mapping between local numbering and global numbering
912: .  nrows - the number of rows indicated
913: .  rows - the indices of the rows
914: .  ncols - the number of columns in the matrix
915: .  cols - the columns indicated
916: .  dnz - the array that will be passed to the matrix preallocation routines
917: -  onz - the other array passed to the matrix preallocation routines

919:    Level: intermediate

921:    Notes:
922:     See Users-Manual: ch_performance for more details.

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

926: .seealso: MatPreallocateFinalize(), MatPreallocateSet()
927:           MatPreallocateInitialize(),  MatPreallocateSetLocal()
928: M*/
929: #define MatPreallocateSymmetricSetLocalBlock(map,nrows,rows,ncols,cols,dnz,onz) 0;\
930: do {\
931:   PetscInt __l;\
932:   _4_ISLocalToGlobalMappingApplyBlock(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
933:   _4_ISLocalToGlobalMappingApplyBlock(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
934:   for (__l=0;__l<nrows;__l++) {\
935:     _4_MatPreallocateSymmetricSetBlock((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
936:   }\
937: } while (0)

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

943:    Synopsis:
944: #include <petscmat.h>
945:    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)

947:    Not Collective

949:    Input Parameters:
950: +  row - the row
951: .  ncols - the number of columns in the matrix
952: -  cols - the columns indicated

954:    Output Parameters:
955: +  dnz - the array that will be passed to the matrix preallocation routines
956: -  onz - the other array passed to the matrix preallocation routines

958:    Level: intermediate

960:    Notes:
961:     See Users-Manual: ch_performance for more details.

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

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

967: .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock()
968:           MatPreallocateInitialize(), MatPreallocateSetLocal()
969: M*/
970: #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
971: do { PetscInt __i; \
972:   if (row < __rstart) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Trying to set preallocation for row %D less than first local row %D",row,__rstart);\
973:   if (row >= __rstart+__nrows) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Trying to set preallocation for row %D greater than last local row %D",row,__rstart+__nrows-1);\
974:   for (__i=0; __i<nc; __i++) {\
975:     if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \
976:     else if (dnz[row - __rstart] < __ncols) dnz[row - __rstart]++;\
977:   }\
978: } while (0)

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

984:    Synopsis:
985: #include <petscmat.h>
986:    PetscErrorCode MatPreallocateSymmetricSetBlock(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)

988:    Not Collective

990:    Input Parameters:
991: +  nrows - the number of rows indicated
992: .  rows - the indices of the rows
993: .  ncols - the number of columns in the matrix
994: .  cols - the columns indicated
995: .  dnz - the array that will be passed to the matrix preallocation routines
996: -  onz - the other array passed to the matrix preallocation routines

998:    Level: intermediate

1000:    Notes:
1001:     See Users-Manual: ch_performance for more details.

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

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

1007: .seealso: MatPreallocateFinalize(), MatPreallocateSet(),  MatPreallocateInitialize(),
1008:           MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal()
1009: M*/
1010: #define MatPreallocateSymmetricSetBlock(row,nc,cols,dnz,onz) 0;\
1011: do { PetscInt __i; \
1012:   for (__i=0; __i<nc; __i++) {\
1013:     if (cols[__i] >= __end) onz[row - __rstart]++; \
1014:     else if (cols[__i] >= row && dnz[row - __rstart] < __ncols) dnz[row - __rstart]++;\
1015:   }\
1016: } while (0)

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

1021:    Synopsis:
1022: #include <petscmat.h>
1023:    PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz)

1025:    Not Collective

1027:    Input Parameters:
1028: +  A - matrix
1029: .  row - row where values exist (must be local to this process)
1030: .  ncols - number of columns
1031: .  cols - columns with nonzeros
1032: .  dnz - the array that will be passed to the matrix preallocation routines
1033: -  onz - the other array passed to the matrix preallocation routines

1035:    Level: intermediate

1037:    Notes:
1038:     See Users-Manual: ch_performance for more details.

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

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

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


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

1054:    Synopsis:
1055: #include <petscmat.h>
1056:    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)

1058:    Collective

1060:    Input Parameters:
1061: +  dnz - the array that was be passed to the matrix preallocation routines
1062: -  onz - the other array passed to the matrix preallocation routines

1064:    Level: intermediate

1066:    Notes:
1067:     See Users-Manual: ch_performance for more details.

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

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

1073: .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
1074:           MatPreallocateSymmetricSetLocalBlock()
1075: M*/
1076: #define MatPreallocateFinalize(dnz,onz) 0;_4_PetscFree2(dnz,onz);CHKERRQ(_4_ierr);} while (0)

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

1081: PETSC_EXTERN PetscErrorCode MatInodeAdjustForInodes(Mat,IS*,IS*);
1082: PETSC_EXTERN PetscErrorCode MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);

1084: PETSC_EXTERN PetscErrorCode MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
1085: PETSC_EXTERN PetscErrorCode MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
1086: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
1087: PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
1088: PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
1089: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJFromTriple(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*,PetscInt,PetscBool);

1091: #define MAT_SKIP_ALLOCATION -4

1093: PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
1094: PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
1095: PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
1096: PETSC_EXTERN PetscErrorCode MatSeqAIJSetTotalPreallocation(Mat,PetscInt);

1098: PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1099: PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1100: PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1101: PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []);
1102: PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
1103: PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
1104: PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
1105: PETSC_EXTERN PetscErrorCode MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
1106: PETSC_EXTERN PetscErrorCode MatMPIAdjToSeq(Mat,Mat*);
1107: PETSC_EXTERN PetscErrorCode MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
1108: PETSC_EXTERN PetscErrorCode MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
1109: PETSC_EXTERN PetscErrorCode MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,const PetscInt*[]);
1110: PETSC_EXTERN PetscErrorCode MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,const PetscInt*[]);
1111: PETSC_EXTERN PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat,Mat*);

1113: PETSC_EXTERN PetscErrorCode MatDenseGetLDA(Mat,PetscInt*);
1114: PETSC_EXTERN PetscErrorCode MatDenseSetLDA(Mat,PetscInt);
1115: PETSC_DEPRECATED_FUNCTION("Use MatDenseSetLDA() (since version 3.14)") PETSC_STATIC_INLINE PetscErrorCode MatSeqDenseSetLDA(Mat A,PetscInt lda) {return MatDenseSetLDA(A,lda);}
1116: PETSC_EXTERN PetscErrorCode MatDenseGetLocalMatrix(Mat,Mat*);

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

1120: PETSC_EXTERN PetscErrorCode MatStoreValues(Mat);
1121: PETSC_EXTERN PetscErrorCode MatRetrieveValues(Mat);

1123: PETSC_EXTERN PetscErrorCode MatFindNonzeroRows(Mat,IS*);
1124: PETSC_EXTERN PetscErrorCode MatFindZeroRows(Mat,IS*);
1125: /*
1126:   These routines are not usually accessed directly, rather solving is
1127:   done through the KSP and PC interfaces.
1128: */

1130: /*J
1131:     MatOrderingType - String with the name of a PETSc matrix ordering

1133:    Level: beginner

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

1138: .seealso: MatGetOrdering()
1139: J*/
1140: typedef const char* MatOrderingType;
1141: #define MATORDERINGNATURAL        "natural"
1142: #define MATORDERINGND             "nd"
1143: #define MATORDERING1WD            "1wd"
1144: #define MATORDERINGRCM            "rcm"
1145: #define MATORDERINGQMD            "qmd"
1146: #define MATORDERINGROWLENGTH      "rowlength"
1147: #define MATORDERINGWBM            "wbm"
1148: #define MATORDERINGSPECTRAL       "spectral"
1149: #define MATORDERINGAMD            "amd"            /* only works if UMFPACK is installed with PETSc */
1150: #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 */
1151: #define MATORDERINGEXTERNAL       "external"       /* uses an ordering type internal to the factorization package */

1153: PETSC_EXTERN PetscErrorCode MatGetOrdering(Mat,MatOrderingType,IS*,IS*);
1154: PETSC_EXTERN PetscErrorCode MatGetOrderingList(PetscFunctionList*);
1155: PETSC_EXTERN PetscErrorCode MatOrderingRegister(const char[],PetscErrorCode(*)(Mat,MatOrderingType,IS*,IS*));
1156: PETSC_EXTERN PetscFunctionList MatOrderingList;

1158: PETSC_EXTERN PetscErrorCode MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
1159: PETSC_EXTERN PetscErrorCode MatCreateLaplacian(Mat,PetscReal,PetscBool,Mat*);

1161: /*S
1162:     MatFactorShiftType - Numeric Shift for factorizations

1164:    Level: beginner

1166: .seealso: MatGetFactor()
1167: S*/
1168: typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType;
1169: PETSC_EXTERN const char *const MatFactorShiftTypes[];
1170: PETSC_EXTERN const char *const MatFactorShiftTypesDetail[];

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

1175:     Level: beginner

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

1180: .seealso: MatGetFactor()
1181: S*/
1182: typedef enum {MAT_FACTOR_NOERROR,MAT_FACTOR_STRUCT_ZEROPIVOT,MAT_FACTOR_NUMERIC_ZEROPIVOT,MAT_FACTOR_OUTMEMORY,MAT_FACTOR_OTHER} MatFactorError;

1184: PETSC_EXTERN PetscErrorCode MatFactorGetError(Mat,MatFactorError*);
1185: PETSC_EXTERN PetscErrorCode MatFactorClearError(Mat);
1186: PETSC_EXTERN PetscErrorCode MatFactorGetErrorZeroPivot(Mat,PetscReal*,PetscInt*);

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

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

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

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

1199:    Level: developer

1201: .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
1202:           MatFactorInfoInitialize()

1204: S*/
1205: typedef struct {
1206:   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
1207:   PetscReal     usedt;
1208:   PetscReal     dt;             /* drop tolerance */
1209:   PetscReal     dtcol;          /* tolerance for pivoting */
1210:   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
1211:   PetscReal     fill;           /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */
1212:   PetscReal     levels;         /* ICC/ILU(levels) */
1213:   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
1214:                                    factorization may be faster if do not pivot */
1215:   PetscReal     zeropivot;      /* pivot is called zero if less than this */
1216:   PetscReal     shifttype;      /* type of shift added to matrix factor to prevent zero pivots */
1217:   PetscReal     shiftamount;     /* how large the shift is */
1218: } MatFactorInfo;

1220: PETSC_EXTERN PetscErrorCode MatFactorInfoInitialize(MatFactorInfo*);
1221: PETSC_EXTERN PetscErrorCode MatCholeskyFactor(Mat,IS,const MatFactorInfo*);
1222: PETSC_EXTERN PetscErrorCode MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
1223: PETSC_EXTERN PetscErrorCode MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*);
1224: PETSC_EXTERN PetscErrorCode MatLUFactor(Mat,IS,IS,const MatFactorInfo*);
1225: PETSC_EXTERN PetscErrorCode MatILUFactor(Mat,IS,IS,const MatFactorInfo*);
1226: PETSC_EXTERN PetscErrorCode MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
1227: PETSC_EXTERN PetscErrorCode MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
1228: PETSC_EXTERN PetscErrorCode MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
1229: PETSC_EXTERN PetscErrorCode MatICCFactor(Mat,IS,const MatFactorInfo*);
1230: PETSC_EXTERN PetscErrorCode MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*);
1231: PETSC_EXTERN PetscErrorCode MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
1232: PETSC_EXTERN PetscErrorCode MatSolve(Mat,Vec,Vec);
1233: PETSC_EXTERN PetscErrorCode MatForwardSolve(Mat,Vec,Vec);
1234: PETSC_EXTERN PetscErrorCode MatBackwardSolve(Mat,Vec,Vec);
1235: PETSC_EXTERN PetscErrorCode MatSolveAdd(Mat,Vec,Vec,Vec);
1236: PETSC_EXTERN PetscErrorCode MatSolveTranspose(Mat,Vec,Vec);
1237: PETSC_EXTERN PetscErrorCode MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
1238: PETSC_EXTERN PetscErrorCode MatSolves(Mat,Vecs,Vecs);
1239: PETSC_EXTERN PetscErrorCode MatSetUnfactored(Mat);

1241: typedef enum {MAT_FACTOR_SCHUR_UNFACTORED, MAT_FACTOR_SCHUR_FACTORED, MAT_FACTOR_SCHUR_INVERTED} MatFactorSchurStatus;
1242: PETSC_EXTERN PetscErrorCode MatFactorSetSchurIS(Mat,IS);
1243: PETSC_EXTERN PetscErrorCode MatFactorGetSchurComplement(Mat,Mat*,MatFactorSchurStatus*);
1244: PETSC_EXTERN PetscErrorCode MatFactorRestoreSchurComplement(Mat,Mat*,MatFactorSchurStatus);
1245: PETSC_EXTERN PetscErrorCode MatFactorInvertSchurComplement(Mat);
1246: PETSC_EXTERN PetscErrorCode MatFactorCreateSchurComplement(Mat,Mat*,MatFactorSchurStatus*);
1247: PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplement(Mat,Vec,Vec);
1248: PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplementTranspose(Mat,Vec,Vec);
1249: PETSC_EXTERN PetscErrorCode MatFactorFactorizeSchurComplement(Mat);

1251: /*E
1252:     MatSORType - What type of (S)SOR to perform

1254:     Level: beginner

1256:    May be bitwise ORd together

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

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

1262: .seealso: MatSOR()
1263: E*/
1264: typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1265:               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1266:               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
1267:               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
1268: PETSC_EXTERN PetscErrorCode MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);

1270: /*
1271:     These routines are for efficiently computing Jacobians via finite differences.
1272: */

1274: /*S
1275:      MatColoring - Object for managing the coloring of matrices.

1277:    Level: beginner

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

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

1286: .seealso:  MatFDColoringCreate(), MatColoringWeightType, ISColoring, MatFDColoring, DMCreateColoring(), MatColoringCreate(), MatOrdering, MatPartitioning
1287: S*/
1288: typedef struct _p_MatColoring* MatColoring;

1290: /*J
1291:     MatColoringType - String with the name of a PETSc matrix coloring

1293:    Level: beginner

1295: .seealso: MatColoringSetType(), MatColoring
1296: J*/
1297: typedef const  char*           MatColoringType;
1298: #define MATCOLORINGJP      "jp"
1299: #define MATCOLORINGPOWER   "power"
1300: #define MATCOLORINGNATURAL "natural"
1301: #define MATCOLORINGSL      "sl"
1302: #define MATCOLORINGLF      "lf"
1303: #define MATCOLORINGID      "id"
1304: #define MATCOLORINGGREEDY  "greedy"

1306: /*E
1307:    MatColoringWeightType - Type of weight scheme

1309:     Not Collective

1311: +   MAT_COLORING_RANDOM  - Random weights
1312: .   MAT_COLORING_LEXICAL - Lexical weighting based upon global numbering.
1313: -   MAT_COLORING_LF      - Last-first weighting.

1315:     Level: intermediate

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

1319: .seealso: MatColoring, MatColoringCreate()
1320: E*/
1321: typedef enum {MAT_COLORING_WEIGHT_RANDOM,MAT_COLORING_WEIGHT_LEXICAL,MAT_COLORING_WEIGHT_LF,MAT_COLORING_WEIGHT_SL} MatColoringWeightType;

1323: PETSC_EXTERN PetscErrorCode MatColoringCreate(Mat,MatColoring*);
1324: PETSC_EXTERN PetscErrorCode MatColoringGetDegrees(Mat,PetscInt,PetscInt*);
1325: PETSC_EXTERN PetscErrorCode MatColoringDestroy(MatColoring*);
1326: PETSC_EXTERN PetscErrorCode MatColoringView(MatColoring,PetscViewer);
1327: PETSC_EXTERN PetscErrorCode MatColoringSetType(MatColoring,MatColoringType);
1328: PETSC_EXTERN PetscErrorCode MatColoringSetFromOptions(MatColoring);
1329: PETSC_EXTERN PetscErrorCode MatColoringSetDistance(MatColoring,PetscInt);
1330: PETSC_EXTERN PetscErrorCode MatColoringGetDistance(MatColoring,PetscInt*);
1331: PETSC_EXTERN PetscErrorCode MatColoringSetMaxColors(MatColoring,PetscInt);
1332: PETSC_EXTERN PetscErrorCode MatColoringGetMaxColors(MatColoring,PetscInt*);
1333: PETSC_EXTERN PetscErrorCode MatColoringApply(MatColoring,ISColoring*);
1334: PETSC_EXTERN PetscErrorCode MatColoringRegister(const char[],PetscErrorCode(*)(MatColoring));
1335: PETSC_EXTERN PetscErrorCode MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
1336: PETSC_EXTERN PetscErrorCode MatColoringSetWeightType(MatColoring,MatColoringWeightType);
1337: PETSC_EXTERN PetscErrorCode MatColoringSetWeights(MatColoring,PetscReal*,PetscInt*);
1338: PETSC_EXTERN PetscErrorCode MatColoringCreateWeights(MatColoring,PetscReal **,PetscInt **lperm);
1339: PETSC_EXTERN PetscErrorCode MatColoringTest(MatColoring,ISColoring);
1340: PETSC_DEPRECATED_FUNCTION("Use MatColoringTest() (since version 3.10)") PETSC_STATIC_INLINE PetscErrorCode MatColoringTestValid(MatColoring matcoloring,ISColoring iscoloring) {return MatColoringTest(matcoloring,iscoloring);}
1341: PETSC_EXTERN PetscErrorCode MatISColoringTest(Mat,ISColoring);

1343: /*S
1344:      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1345:         and coloring

1347:    Level: beginner

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

1352: .seealso:  MatFDColoringCreate(), MatColoring, DMCreateColoring()
1353: S*/
1354: typedef struct _p_MatFDColoring* MatFDColoring;

1356: PETSC_EXTERN PetscErrorCode MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
1357: PETSC_EXTERN PetscErrorCode MatFDColoringDestroy(MatFDColoring*);
1358: PETSC_EXTERN PetscErrorCode MatFDColoringView(MatFDColoring,PetscViewer);
1359: PETSC_EXTERN PetscErrorCode MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
1360: PETSC_EXTERN PetscErrorCode MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**);
1361: PETSC_EXTERN PetscErrorCode MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
1362: PETSC_EXTERN PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring);
1363: PETSC_EXTERN PetscErrorCode MatFDColoringApply(Mat,MatFDColoring,Vec,void *);
1364: PETSC_EXTERN PetscErrorCode MatFDColoringSetF(MatFDColoring,Vec);
1365: PETSC_EXTERN PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,const PetscInt*[]);
1366: PETSC_EXTERN PetscErrorCode MatFDColoringSetUp(Mat,ISColoring,MatFDColoring);
1367: PETSC_EXTERN PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring,PetscInt,PetscInt);
1368: PETSC_EXTERN PetscErrorCode MatFDColoringSetValues(Mat,MatFDColoring,const PetscScalar*);

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

1373:    Level: beginner

1375: .seealso:  MatTransposeColoringCreate()
1376: S*/
1377: typedef struct _p_MatTransposeColoring* MatTransposeColoring;

1379: PETSC_EXTERN PetscErrorCode MatTransposeColoringCreate(Mat,ISColoring,MatTransposeColoring *);
1380: PETSC_EXTERN PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring,Mat,Mat);
1381: PETSC_EXTERN PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring,Mat,Mat);
1382: PETSC_EXTERN PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring*);

1384: /*
1385:     These routines are for partitioning matrices: currently used only
1386:   for adjacency matrix, MatCreateMPIAdj().
1387: */

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

1392:    Level: beginner

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

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

1401: .seealso:  MatPartitioningCreate(), MatPartitioningType, MatColoring, MatGetOrdering()
1402: S*/
1403: typedef struct _p_MatPartitioning* MatPartitioning;

1405: /*J
1406:     MatPartitioningType - String with the name of a PETSc matrix partitioning

1408:    Level: beginner
1409: dm
1410: .seealso: MatPartitioningCreate(), MatPartitioning
1411: J*/
1412: typedef const char* MatPartitioningType;
1413: #define MATPARTITIONINGCURRENT  "current"
1414: #define MATPARTITIONINGAVERAGE   "average"
1415: #define MATPARTITIONINGSQUARE   "square"
1416: #define MATPARTITIONINGPARMETIS "parmetis"
1417: #define MATPARTITIONINGCHACO    "chaco"
1418: #define MATPARTITIONINGPARTY    "party"
1419: #define MATPARTITIONINGPTSCOTCH "ptscotch"
1420: #define MATPARTITIONINGHIERARCH  "hierarch"

1422: PETSC_EXTERN PetscErrorCode MatPartitioningCreate(MPI_Comm,MatPartitioning*);
1423: PETSC_EXTERN PetscErrorCode MatPartitioningSetType(MatPartitioning,MatPartitioningType);
1424: PETSC_EXTERN PetscErrorCode MatPartitioningSetNParts(MatPartitioning,PetscInt);
1425: PETSC_EXTERN PetscErrorCode MatPartitioningSetAdjacency(MatPartitioning,Mat);
1426: PETSC_EXTERN PetscErrorCode MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
1427: PETSC_EXTERN PetscErrorCode MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
1428: PETSC_EXTERN PetscErrorCode MatPartitioningSetUseEdgeWeights(MatPartitioning,PetscBool);
1429: PETSC_EXTERN PetscErrorCode MatPartitioningGetUseEdgeWeights(MatPartitioning,PetscBool*);
1430: PETSC_EXTERN PetscErrorCode MatPartitioningApply(MatPartitioning,IS*);
1431: PETSC_EXTERN PetscErrorCode MatPartitioningImprove(MatPartitioning,IS*);
1432: PETSC_EXTERN PetscErrorCode MatPartitioningViewImbalance(MatPartitioning,IS);
1433: PETSC_EXTERN PetscErrorCode MatPartitioningApplyND(MatPartitioning,IS*);
1434: PETSC_EXTERN PetscErrorCode MatPartitioningDestroy(MatPartitioning*);
1435: PETSC_EXTERN PetscErrorCode MatPartitioningRegister(const char[],PetscErrorCode (*)(MatPartitioning));
1436: PETSC_EXTERN PetscErrorCode MatPartitioningView(MatPartitioning,PetscViewer);
1437: PETSC_EXTERN PetscErrorCode MatPartitioningViewFromOptions(MatPartitioning,PetscObject,const char[]);
1438: PETSC_EXTERN PetscErrorCode MatPartitioningSetFromOptions(MatPartitioning);
1439: PETSC_EXTERN PetscErrorCode MatPartitioningGetType(MatPartitioning,MatPartitioningType*);

1441: PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetRepartition(MatPartitioning part);
1442: PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
1443: PETSC_EXTERN PetscErrorCode MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *);

1445: typedef enum { MP_CHACO_MULTILEVEL=1,MP_CHACO_SPECTRAL=2,MP_CHACO_LINEAR=4,MP_CHACO_RANDOM=5,MP_CHACO_SCATTERED=6 } MPChacoGlobalType;
1446: PETSC_EXTERN const char *const MPChacoGlobalTypes[];
1447: typedef enum { MP_CHACO_KERNIGHAN=1,MP_CHACO_NONE=2 } MPChacoLocalType;
1448: PETSC_EXTERN const char *const MPChacoLocalTypes[];
1449: typedef enum { MP_CHACO_LANCZOS=0,MP_CHACO_RQI=1 } MPChacoEigenType;
1450: PETSC_EXTERN const char *const MPChacoEigenTypes[];

1452: PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetGlobal(MatPartitioning,MPChacoGlobalType);
1453: PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetGlobal(MatPartitioning,MPChacoGlobalType*);
1454: PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetLocal(MatPartitioning,MPChacoLocalType);
1455: PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetLocal(MatPartitioning,MPChacoLocalType*);
1456: PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
1457: PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
1458: PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenSolver(MatPartitioning,MPChacoEigenType*);
1459: PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenTol(MatPartitioning,PetscReal);
1460: PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenTol(MatPartitioning,PetscReal*);
1461: PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenNumber(MatPartitioning,PetscInt);
1462: PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenNumber(MatPartitioning,PetscInt*);

1464: #define MP_PARTY_OPT "opt"
1465: #define MP_PARTY_LIN "lin"
1466: #define MP_PARTY_SCA "sca"
1467: #define MP_PARTY_RAN "ran"
1468: #define MP_PARTY_GBF "gbf"
1469: #define MP_PARTY_GCF "gcf"
1470: #define MP_PARTY_BUB "bub"
1471: #define MP_PARTY_DEF "def"
1472: PETSC_EXTERN PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning,const char*);
1473: #define MP_PARTY_HELPFUL_SETS "hs"
1474: #define MP_PARTY_KERNIGHAN_LIN "kl"
1475: #define MP_PARTY_NONE "no"
1476: PETSC_EXTERN PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning,const char*);
1477: PETSC_EXTERN PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
1478: PETSC_EXTERN PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning,PetscBool);
1479: PETSC_EXTERN PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool);

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

1484: PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning,PetscReal);
1485: PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning,PetscReal*);
1486: PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning,MPPTScotchStrategyType);
1487: PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning,MPPTScotchStrategyType*);

1489: /*
1490:  * hierarchical partitioning
1491:  */
1492: PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetFineparts(MatPartitioning,IS*);
1493: PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetCoarseparts(MatPartitioning,IS*);
1494: PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNcoarseparts(MatPartitioning,PetscInt);
1495: PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNfineparts(MatPartitioning, PetscInt);

1497: PETSC_EXTERN PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*);
1498: PETSC_EXTERN PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*);

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

1667: /*
1668:    Codes for matrices stored on disk. By default they are
1669:    stored in a universal format. By changing the format with
1670:    PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will
1671:    be stored in a way natural for the matrix, for example dense matrices
1672:    would be stored as dense. Matrices stored this way may only be
1673:    read into matrices of the same type.
1674: */
1675: #define MATRIX_BINARY_FORMAT_DENSE -1

1677: PETSC_EXTERN PetscErrorCode MatMPIBAIJSetHashTableFactor(Mat,PetscReal);

1679: PETSC_EXTERN PetscErrorCode MatISSetLocalMatType(Mat,MatType);
1680: PETSC_EXTERN PetscErrorCode MatISSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1681: PETSC_EXTERN PetscErrorCode MatISStoreL2L(Mat,PetscBool);
1682: PETSC_EXTERN PetscErrorCode MatISFixLocalEmpty(Mat,PetscBool);
1683: PETSC_EXTERN PetscErrorCode MatISGetLocalMat(Mat,Mat*);
1684: PETSC_EXTERN PetscErrorCode MatISRestoreLocalMat(Mat,Mat*);
1685: PETSC_EXTERN PetscErrorCode MatISSetLocalMat(Mat,Mat);
1686: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use the MatConvert() interface (since version 3.10)") PetscErrorCode MatISGetMPIXAIJ(Mat,MatReuse,Mat*);

1688: /*S
1689:      MatNullSpace - Object that removes a null space from a vector, i.e.
1690:          orthogonalizes the vector to a subsapce

1692:    Level: advanced

1694:   Users manual sections:
1695: .   sec_singular

1697: .seealso:  MatNullSpaceCreate()
1698: S*/
1699: typedef struct _p_MatNullSpace* MatNullSpace;

1701: PETSC_EXTERN PetscErrorCode MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*);
1702: PETSC_EXTERN PetscErrorCode MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*);
1703: PETSC_EXTERN PetscErrorCode MatNullSpaceDestroy(MatNullSpace*);
1704: PETSC_EXTERN PetscErrorCode MatNullSpaceRemove(MatNullSpace,Vec);
1705: PETSC_EXTERN PetscErrorCode MatGetNullSpace(Mat, MatNullSpace *);
1706: PETSC_EXTERN PetscErrorCode MatGetTransposeNullSpace(Mat, MatNullSpace *);
1707: PETSC_EXTERN PetscErrorCode MatSetTransposeNullSpace(Mat,MatNullSpace);
1708: PETSC_EXTERN PetscErrorCode MatSetNullSpace(Mat,MatNullSpace);
1709: PETSC_EXTERN PetscErrorCode MatSetNearNullSpace(Mat,MatNullSpace);
1710: PETSC_EXTERN PetscErrorCode MatGetNearNullSpace(Mat,MatNullSpace*);
1711: PETSC_EXTERN PetscErrorCode MatNullSpaceTest(MatNullSpace,Mat,PetscBool  *);
1712: PETSC_EXTERN PetscErrorCode MatNullSpaceView(MatNullSpace,PetscViewer);
1713: PETSC_EXTERN PetscErrorCode MatNullSpaceGetVecs(MatNullSpace,PetscBool*,PetscInt*,const Vec**);
1714: PETSC_EXTERN PetscErrorCode MatNullSpaceCreateRigidBody(Vec,MatNullSpace*);

1716: PETSC_EXTERN PetscErrorCode MatReorderingSeqSBAIJ(Mat,IS);
1717: PETSC_EXTERN PetscErrorCode MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
1718: PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);

1720: PETSC_EXTERN PetscErrorCode MatCreateMAIJ(Mat,PetscInt,Mat*);
1721: PETSC_EXTERN PetscErrorCode MatMAIJRedimension(Mat,PetscInt,Mat*);
1722: PETSC_EXTERN PetscErrorCode MatMAIJGetAIJ(Mat,Mat*);

1724: PETSC_EXTERN PetscErrorCode MatComputeOperator(Mat,MatType,Mat*);
1725: PETSC_EXTERN PetscErrorCode MatComputeOperatorTranspose(Mat,MatType,Mat*);

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

1730: PETSC_EXTERN PetscErrorCode MatCreateKAIJ(Mat,PetscInt,PetscInt,const PetscScalar[],const PetscScalar[],Mat*);
1731: PETSC_EXTERN PetscErrorCode MatKAIJGetAIJ(Mat,Mat*);
1732: PETSC_EXTERN PetscErrorCode MatKAIJGetS(Mat,PetscInt*,PetscInt*,PetscScalar**);
1733: PETSC_EXTERN PetscErrorCode MatKAIJGetSRead(Mat,PetscInt*,PetscInt*,const PetscScalar**);
1734: PETSC_EXTERN PetscErrorCode MatKAIJRestoreS(Mat,PetscScalar**);
1735: PETSC_EXTERN PetscErrorCode MatKAIJRestoreSRead(Mat,const PetscScalar**);
1736: PETSC_EXTERN PetscErrorCode MatKAIJGetT(Mat,PetscInt*,PetscInt*,PetscScalar**);
1737: PETSC_EXTERN PetscErrorCode MatKAIJGetTRead(Mat,PetscInt*,PetscInt*,const PetscScalar**);
1738: PETSC_EXTERN PetscErrorCode MatKAIJRestoreT(Mat,PetscScalar**);
1739: PETSC_EXTERN PetscErrorCode MatKAIJRestoreTRead(Mat,const PetscScalar**);
1740: PETSC_EXTERN PetscErrorCode MatKAIJSetAIJ(Mat,Mat);
1741: PETSC_EXTERN PetscErrorCode MatKAIJSetS(Mat,PetscInt,PetscInt,const PetscScalar[]);
1742: PETSC_EXTERN PetscErrorCode MatKAIJSetT(Mat,PetscInt,PetscInt,const PetscScalar[]);
1743: PETSC_EXTERN PetscErrorCode MatKAIJGetScaledIdentity(Mat,PetscBool*);

1745: PETSC_EXTERN PetscErrorCode MatDiagonalScaleLocal(Mat,Vec);

1747: PETSC_EXTERN PetscErrorCode MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
1748: PETSC_EXTERN PetscErrorCode MatMFFDSetBase(Mat,Vec,Vec);
1749: PETSC_EXTERN PetscErrorCode MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*);
1750: PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*));
1751: PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec));
1752: PETSC_EXTERN PetscErrorCode MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt);
1753: PETSC_EXTERN PetscErrorCode MatMFFDResetHHistory(Mat);
1754: PETSC_EXTERN PetscErrorCode MatMFFDSetFunctionError(Mat,PetscReal);
1755: PETSC_EXTERN PetscErrorCode MatMFFDSetPeriod(Mat,PetscInt);
1756: PETSC_EXTERN PetscErrorCode MatMFFDGetH(Mat,PetscScalar *);
1757: PETSC_EXTERN PetscErrorCode MatMFFDSetOptionsPrefix(Mat,const char[]);
1758: PETSC_EXTERN PetscErrorCode MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*);
1759: PETSC_EXTERN PetscErrorCode MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*);

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

1765:     Notes:
1766:     MATMFFD is a specific MatType which uses the MatMFFD data structure

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

1770:     Level: developer

1772: .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister()
1773: S*/
1774: typedef struct _p_MatMFFD* MatMFFD;

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

1779:    Level: beginner

1781: .seealso: MatMFFDSetType(), MatMFFDRegister()
1782: J*/
1783: typedef const char* MatMFFDType;
1784: #define MATMFFD_DS  "ds"
1785: #define MATMFFD_WP  "wp"

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

1790: PETSC_EXTERN PetscErrorCode MatMFFDDSSetUmin(Mat,PetscReal);
1791: PETSC_EXTERN PetscErrorCode MatMFFDWPSetComputeNormU(Mat,PetscBool);

1793: PETSC_EXTERN PetscErrorCode MatFDColoringSetType(MatFDColoring,MatMFFDType);

1795: PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
1796: PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);

1798: #ifdef PETSC_HAVE_HARA
1799: PETSC_EXTERN_TYPEDEF typedef PetscScalar (*MatHaraKernel)(PetscInt,PetscReal[],PetscReal[],void*);
1800: PETSC_EXTERN PetscErrorCode MatCreateHaraFromKernel(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscReal[],MatHaraKernel,void*,PetscReal,PetscInt,PetscInt,Mat*);
1801: PETSC_EXTERN PetscErrorCode MatCreateHaraFromMat(Mat,PetscInt,const PetscReal[],PetscReal,PetscInt,PetscInt,PetscInt,PetscReal,Mat*);
1802: PETSC_EXTERN PetscErrorCode MatHaraSetSamplingMat(Mat,Mat,PetscInt,PetscReal);
1803: #endif

1805: /*
1806:    PETSc interface to MUMPS
1807: */
1808: #ifdef PETSC_HAVE_MUMPS
1809: PETSC_EXTERN PetscErrorCode MatMumpsSetIcntl(Mat,PetscInt,PetscInt);
1810: PETSC_EXTERN PetscErrorCode MatMumpsGetIcntl(Mat,PetscInt,PetscInt*);
1811: PETSC_EXTERN PetscErrorCode MatMumpsSetCntl(Mat,PetscInt,PetscReal);
1812: PETSC_EXTERN PetscErrorCode MatMumpsGetCntl(Mat,PetscInt,PetscReal*);

1814: PETSC_EXTERN PetscErrorCode MatMumpsGetInfo(Mat,PetscInt,PetscInt*);
1815: PETSC_EXTERN PetscErrorCode MatMumpsGetInfog(Mat,PetscInt,PetscInt*);
1816: PETSC_EXTERN PetscErrorCode MatMumpsGetRinfo(Mat,PetscInt,PetscReal*);
1817: PETSC_EXTERN PetscErrorCode MatMumpsGetRinfog(Mat,PetscInt,PetscReal*);
1818: PETSC_EXTERN PetscErrorCode MatMumpsGetInverse(Mat,Mat);
1819: PETSC_EXTERN PetscErrorCode MatMumpsGetInverseTranspose(Mat,Mat);
1820: #endif

1822: /*
1823:    PETSc interface to Mkl_Pardiso
1824: */
1825: #ifdef PETSC_HAVE_MKL_PARDISO
1826: PETSC_EXTERN PetscErrorCode MatMkl_PardisoSetCntl(Mat,PetscInt,PetscInt);
1827: #endif

1829: /*
1830:    PETSc interface to Mkl_CPardiso
1831: */
1832: #ifdef PETSC_HAVE_MKL_CPARDISO
1833: PETSC_EXTERN PetscErrorCode MatMkl_CPardisoSetCntl(Mat,PetscInt,PetscInt);
1834: #endif

1836: /*
1837:    PETSc interface to SUPERLU
1838: */
1839: #ifdef PETSC_HAVE_SUPERLU
1840: PETSC_EXTERN PetscErrorCode MatSuperluSetILUDropTol(Mat,PetscReal);
1841: #endif

1843: /*
1844:    PETSc interface to SUPERLU_DIST
1845: */
1846: #ifdef PETSC_HAVE_SUPERLU_DIST
1847: PETSC_EXTERN PetscErrorCode MatSuperluDistGetDiagU(Mat,PetscScalar*);
1848: #endif

1850: /*
1851:    PETSc interface to STRUMPACK
1852: */
1853: #ifdef PETSC_HAVE_STRUMPACK
1854: /*E
1855:     MatSTRUMPACKReordering - sparsity reducing ordering to be used in STRUMPACK

1857:     Level: intermediate
1858: E*/
1859: typedef enum {MAT_STRUMPACK_NATURAL=0,
1860:               MAT_STRUMPACK_METIS=1,
1861:               MAT_STRUMPACK_PARMETIS=2,
1862:               MAT_STRUMPACK_SCOTCH=3,
1863:               MAT_STRUMPACK_PTSCOTCH=4,
1864:               MAT_STRUMPACK_RCM=5} MatSTRUMPACKReordering;
1865: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetReordering(Mat,MatSTRUMPACKReordering);
1866: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetColPerm(Mat,PetscBool);
1867: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSRelTol(Mat,PetscReal);
1868: PETSC_DEPRECATED_FUNCTION("Use MatSTRUMPACKSetHSSRelTol() (since version 3.9)") PETSC_STATIC_INLINE PetscErrorCode MatSTRUMPACKSetHSSRelCompTol(Mat mat,PetscReal rtol) {return MatSTRUMPACKSetHSSRelTol(mat,rtol);}
1869: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSAbsTol(Mat,PetscReal);
1870: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSMinSepSize(Mat,PetscInt);
1871: PETSC_DEPRECATED_FUNCTION("Use MatSTRUMPACKSetHSSMinSepSize() (since version 3.9)") PETSC_STATIC_INLINE PetscErrorCode MatSTRUMPACKSetHSSMinSize(Mat mat,PetscInt hssminsize) {return MatSTRUMPACKSetHSSMinSepSize(mat,hssminsize);}
1872: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSMaxRank(Mat,PetscInt);
1873: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSLeafSize(Mat,PetscInt);
1874: #endif


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

1880: #ifdef PETSC_HAVE_CUDA
1881: /*E
1882:     MatCUSPARSEStorageFormat - indicates the storage format for CUSPARSE (GPU)
1883:     matrices.

1885:     Not Collective

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

1891:     Level: intermediate

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

1895: .seealso: MatCUSPARSESetFormat(), MatCUSPARSEFormatOperation
1896: E*/

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

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

1903: /*E
1904:     MatCUSPARSEFormatOperation - indicates the operation of CUSPARSE (GPU)
1905:     matrices whose operation should use a particular storage format.

1907:     Not Collective

1909: +   MAT_CUSPARSE_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult
1910: .   MAT_CUSPARSE_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult
1911: .   MAT_CUSPARSE_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult
1912: -   MAT_CUSPARSE_ALL - sets the storage format for all CUSPARSE (GPU) matrices

1914:     Level: intermediate

1916: .seealso: MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat
1917: E*/
1918: typedef enum {MAT_CUSPARSE_MULT_DIAG, MAT_CUSPARSE_MULT_OFFDIAG, MAT_CUSPARSE_MULT, MAT_CUSPARSE_ALL} MatCUSPARSEFormatOperation;

1920: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
1921: PETSC_EXTERN PetscErrorCode MatCreateAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1922: PETSC_EXTERN PetscErrorCode MatCUSPARSESetFormat(Mat,MatCUSPARSEFormatOperation,MatCUSPARSEStorageFormat);
1923: PETSC_EXTERN PetscErrorCode MatSeqAIJCUSPARSESetGenerateTranspose(Mat,PetscBool);

1925: typedef struct _p_SplitCSRMat PetscSplitCSRDataStructure;
1926: PETSC_EXTERN PetscErrorCode MatCUSPARSEGetDeviceMatWrite(Mat,PetscSplitCSRDataStructure**);

1928: PETSC_EXTERN PetscErrorCode MatCreateDenseCUDA(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
1929: PETSC_EXTERN PetscErrorCode MatCreateSeqDenseCUDA(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
1930: PETSC_EXTERN PetscErrorCode MatMPIDenseCUDASetPreallocation(Mat,PetscScalar[]);
1931: PETSC_EXTERN PetscErrorCode MatSeqDenseCUDASetPreallocation(Mat,PetscScalar[]);
1932: PETSC_EXTERN PetscErrorCode MatDenseCUDAGetArrayWrite(Mat,PetscScalar**);
1933: PETSC_EXTERN PetscErrorCode MatDenseCUDAGetArrayRead(Mat,const PetscScalar**);
1934: PETSC_EXTERN PetscErrorCode MatDenseCUDAGetArray(Mat,PetscScalar**);
1935: PETSC_EXTERN PetscErrorCode MatDenseCUDARestoreArrayWrite(Mat,PetscScalar**);
1936: PETSC_EXTERN PetscErrorCode MatDenseCUDARestoreArrayRead(Mat,const PetscScalar**);
1937: PETSC_EXTERN PetscErrorCode MatDenseCUDARestoreArray(Mat,PetscScalar**);
1938: PETSC_EXTERN PetscErrorCode MatDenseCUDAPlaceArray(Mat,const PetscScalar*);
1939: PETSC_EXTERN PetscErrorCode MatDenseCUDAReplaceArray(Mat,const PetscScalar*);
1940: PETSC_EXTERN PetscErrorCode MatDenseCUDAResetArray(Mat);

1942: #endif

1944: #if defined(PETSC_HAVE_VIENNACL)
1945: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
1946: PETSC_EXTERN PetscErrorCode MatCreateAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1947: #endif

1949: /*
1950:    PETSc interface to FFTW
1951: */
1952: #if defined(PETSC_HAVE_FFTW)
1953: PETSC_EXTERN PetscErrorCode VecScatterPetscToFFTW(Mat,Vec,Vec);
1954: PETSC_EXTERN PetscErrorCode VecScatterFFTWToPetsc(Mat,Vec,Vec);
1955: PETSC_EXTERN PetscErrorCode MatCreateVecsFFTW(Mat,Vec*,Vec*,Vec*);
1956: #endif

1958: #if defined(PETSC_HAVE_SCALAPACK)
1959: PETSC_EXTERN PetscErrorCode MatCreateScaLAPACK(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
1960: PETSC_EXTERN PetscErrorCode MatScaLAPACKSetBlockSizes(Mat,PetscInt,PetscInt);
1961: PETSC_EXTERN PetscErrorCode MatScaLAPACKGetBlockSizes(Mat,PetscInt*,PetscInt*);
1962: #endif

1964: PETSC_EXTERN PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*);
1965: PETSC_EXTERN PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*);
1966: PETSC_EXTERN PetscErrorCode MatNestGetISs(Mat,IS[],IS[]);
1967: PETSC_EXTERN PetscErrorCode MatNestGetLocalISs(Mat,IS[],IS[]);
1968: PETSC_EXTERN PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***);
1969: PETSC_EXTERN PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*);
1970: PETSC_EXTERN PetscErrorCode MatNestSetVecType(Mat,VecType);
1971: PETSC_EXTERN PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]);
1972: PETSC_EXTERN PetscErrorCode MatNestSetSubMat(Mat,PetscInt,PetscInt,Mat);

1974: PETSC_EXTERN PetscErrorCode MatChop(Mat,PetscReal);
1975: PETSC_EXTERN PetscErrorCode MatComputeBandwidth(Mat,PetscReal,PetscInt*);

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

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

1981: PETSC_INTERN PetscErrorCode MatHeaderMerge(Mat,Mat*);
1982: PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat,Mat*);
1983: #endif