Actual source code: petscmat.h

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

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

 12:    Level: beginner

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

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

 21:    Level: beginner

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

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

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

112:    Level: beginner

114: .seealso: MatGetFactor(), PCFactorSetMatSolverType(), PCFactorGetMatSolverType()
115: J*/
116: typedef const char* MatSolverType;
117: #define MATSOLVERSUPERLU          "superlu"
118: #define MATSOLVERSUPERLU_DIST     "superlu_dist"
119: #define MATSOLVERSTRUMPACK        "strumpack"
120: #define MATSOLVERUMFPACK          "umfpack"
121: #define MATSOLVERCHOLMOD          "cholmod"
122: #define MATSOLVERKLU              "klu"
123: #define MATSOLVERSPARSEELEMENTAL  "sparseelemental"
124: #define MATSOLVERELEMENTAL        "elemental"
125: #define MATSOLVERESSL             "essl"
126: #define MATSOLVERLUSOL            "lusol"
127: #define MATSOLVERMUMPS            "mumps"
128: #define MATSOLVERMKL_PARDISO      "mkl_pardiso"
129: #define MATSOLVERMKL_CPARDISO     "mkl_cpardiso"
130: #define MATSOLVERPASTIX           "pastix"
131: #define MATSOLVERMATLAB           "matlab"
132: #define MATSOLVERPETSC            "petsc"
133: #define MATSOLVERBAS              "bas"
134: #define MATSOLVERCUSPARSE         "cusparse"
135: #define MATSOLVERCUDA             "cuda"

137: /*E
138:     MatFactorType - indicates what type of factorization is requested

140:     Level: beginner

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

144: .seealso: MatSolverType, MatGetFactor()
145: E*/
146: typedef enum {MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC,MAT_FACTOR_ILUDT} MatFactorType;
147: PETSC_EXTERN const char *const MatFactorTypes[];

149: PETSC_EXTERN PetscErrorCode MatGetFactor(Mat,MatSolverType,MatFactorType,Mat*);
150: PETSC_EXTERN PetscErrorCode MatGetFactorAvailable(Mat,MatSolverType,MatFactorType,PetscBool *);
151: PETSC_EXTERN PetscErrorCode MatFactorGetSolverType(Mat,MatSolverType*);
152: PETSC_EXTERN PetscErrorCode MatGetFactorType(Mat,MatFactorType*);
153: PETSC_EXTERN PetscErrorCode MatSetFactorType(Mat,MatFactorType);
154: PETSC_EXTERN PetscErrorCode MatSolverTypeRegister(MatSolverType,MatType,MatFactorType,PetscErrorCode(*)(Mat,MatFactorType,Mat*));
155: PETSC_EXTERN PetscErrorCode MatSolverTypeGet(MatSolverType,MatType,MatFactorType,PetscBool*,PetscBool*,PetscErrorCode (**)(Mat,MatFactorType,Mat*));
156: typedef MatSolverType MatSolverPackage PETSC_DEPRECATED_TYPEDEF("Use MatSolverType (since version 3.9)");
157: PETSC_DEPRECATED_FUNCTION("Use MatSolverTypeRegister() (since version 3.9)") PETSC_STATIC_INLINE PetscErrorCode MatSolverPackageRegister(MatSolverType stype,MatType mtype,MatFactorType ftype,PetscErrorCode(*f)(Mat,MatFactorType,Mat*))
158: { return MatSolverTypeRegister(stype,mtype,ftype,f); }
159: PETSC_DEPRECATED_FUNCTION("Use MatSolverTypeGet() (since version 3.9)") PETSC_STATIC_INLINE PetscErrorCode MatSolverPackageGet(MatSolverType stype,MatType mtype,MatFactorType ftype,PetscBool *foundmtype,PetscBool *foundstype,PetscErrorCode(**f)(Mat,MatFactorType,Mat*))
160: { return MatSolverTypeGet(stype,mtype,ftype,foundmtype,foundstype,f); }

162: /* Logging support */
163: #define    MAT_FILE_CLASSID 1211216    /* used to indicate matrices in binary files */
164: PETSC_EXTERN PetscClassId MAT_CLASSID;
165: PETSC_EXTERN PetscClassId MAT_COLORING_CLASSID;
166: PETSC_EXTERN PetscClassId MAT_FDCOLORING_CLASSID;
167: PETSC_EXTERN PetscClassId MAT_TRANSPOSECOLORING_CLASSID;
168: PETSC_EXTERN PetscClassId MAT_PARTITIONING_CLASSID;
169: PETSC_EXTERN PetscClassId MAT_COARSEN_CLASSID;
170: PETSC_EXTERN PetscClassId MAT_NULLSPACE_CLASSID;
171: PETSC_EXTERN PetscClassId MATMFFD_CLASSID;

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

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

182:     Level: beginner

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

186: .seealso: MatCreateSubMatrices(), MatCreateSubMatrix(), MatDestroyMatrices(), MatConvert()
187: E*/
188: typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX,MAT_IGNORE_MATRIX,MAT_INPLACE_MATRIX} MatReuse;

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

194:     Level: beginner

196: .seealso: MatGetSeqNonzerostructure()
197: E*/
198: typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatCreateSubMatrixOption;

200: PETSC_EXTERN PetscErrorCode MatInitializePackage(void);

202: PETSC_EXTERN PetscErrorCode MatCreate(MPI_Comm,Mat*);
203: PETSC_EXTERN PetscErrorCode MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt);
204: PETSC_EXTERN PetscErrorCode MatSetType(Mat,MatType);
205: PETSC_EXTERN PetscErrorCode MatSetFromOptions(Mat);
206: PETSC_STATIC_INLINE PetscErrorCode MatViewFromOptions(Mat A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}
207: PETSC_EXTERN PetscErrorCode MatRegister(const char[],PetscErrorCode(*)(Mat));
208: PETSC_EXTERN PetscErrorCode MatRegisterRootName(const char[],const char[],const char[]);
209: PETSC_EXTERN PetscErrorCode MatSetOptionsPrefix(Mat,const char[]);
210: PETSC_EXTERN PetscErrorCode MatAppendOptionsPrefix(Mat,const char[]);
211: PETSC_EXTERN PetscErrorCode MatGetOptionsPrefix(Mat,const char*[]);
212: PETSC_EXTERN PetscErrorCode MatSetErrorIfFailure(Mat,PetscBool);

214: PETSC_EXTERN PetscFunctionList MatList;
215: PETSC_EXTERN PetscFunctionList MatColoringList;
216: PETSC_EXTERN PetscFunctionList MatPartitioningList;

218: /*E
219:     MatStructure - Indicates if two matrices have the same nonzero structure

221:     Level: beginner

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

225: .seealso: MatCopy(), MatAXPY()
226: E*/
227: typedef enum {DIFFERENT_NONZERO_PATTERN,SUBSET_NONZERO_PATTERN,SAME_NONZERO_PATTERN} MatStructure;

229: #if defined PETSC_HAVE_MKL_SPARSE
230: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJMKL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
231: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJMKL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
232: PETSC_EXTERN PetscErrorCode MatCreateBAIJMKL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
233: PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJMKL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
234: #endif

236: PETSC_EXTERN PetscErrorCode MatCreateSeqSELL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
237: PETSC_EXTERN PetscErrorCode MatCreateSELL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
238: PETSC_EXTERN PetscErrorCode MatSeqSELLSetPreallocation(Mat,PetscInt,const PetscInt[]);
239: PETSC_EXTERN PetscErrorCode MatMPISELLSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);

241: PETSC_EXTERN PetscErrorCode MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
242: PETSC_EXTERN PetscErrorCode MatCreateDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
243: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
244: PETSC_EXTERN PetscErrorCode MatCreateAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
245: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
246: PETSC_EXTERN PetscErrorCode MatUpdateMPIAIJWithArrays(Mat,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
247: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*);
248: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithSeqAIJ(MPI_Comm,Mat,Mat,const PetscInt[],Mat*);

250: PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
251: PETSC_EXTERN PetscErrorCode MatCreateBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
252: PETSC_EXTERN PetscErrorCode MatCreateMPIBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat*);

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

257: PETSC_EXTERN PetscErrorCode MatCreateSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
258: PETSC_EXTERN PetscErrorCode MatCreateMPISBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
259: PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
260: PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
261: PETSC_EXTERN PetscErrorCode MatXAIJSetPreallocation(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscInt[],const PetscInt[]);

263: PETSC_EXTERN PetscErrorCode MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*);
264: PETSC_EXTERN PetscErrorCode MatCreateNormal(Mat,Mat*);
265: PETSC_EXTERN PetscErrorCode MatCreateNormalHermitian(Mat,Mat*);
266: PETSC_EXTERN PetscErrorCode MatCreateLRC(Mat,Mat,Vec,Mat,Mat*);
267: PETSC_EXTERN PetscErrorCode MatLRCGetMats(Mat,Mat*,Mat*,Vec*,Mat*);
268: PETSC_EXTERN PetscErrorCode MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,ISLocalToGlobalMapping,Mat*);
269: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
270: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);

272: PETSC_EXTERN PetscErrorCode MatCreateScatter(MPI_Comm,VecScatter,Mat*);
273: PETSC_EXTERN PetscErrorCode MatScatterSetVecScatter(Mat,VecScatter);
274: PETSC_EXTERN PetscErrorCode MatScatterGetVecScatter(Mat,VecScatter*);
275: PETSC_EXTERN PetscErrorCode MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*);
276: PETSC_EXTERN PetscErrorCode MatCompositeAddMat(Mat,Mat);
277: PETSC_EXTERN PetscErrorCode MatCompositeMerge(Mat);
278: typedef enum {MAT_COMPOSITE_MERGE_RIGHT,MAT_COMPOSITE_MERGE_LEFT} MatCompositeMergeType;
279: PETSC_EXTERN PetscErrorCode MatCompositeSetMergeType(Mat,MatCompositeMergeType);
280: PETSC_EXTERN PetscErrorCode MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*);
281: typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType;
282: PETSC_EXTERN PetscErrorCode MatCompositeSetType(Mat,MatCompositeType);
283: PETSC_EXTERN PetscErrorCode MatCompositeGetType(Mat,MatCompositeType*);
284: PETSC_EXTERN PetscErrorCode MatCompositeSetMatStructure(Mat,MatStructure);
285: PETSC_EXTERN PetscErrorCode MatCompositeGetMatStructure(Mat,MatStructure*);
286: PETSC_EXTERN PetscErrorCode MatCompositeGetNumberMat(Mat,PetscInt*);
287: PETSC_EXTERN PetscErrorCode MatCompositeGetMat(Mat,PetscInt,Mat*);
288: PETSC_EXTERN PetscErrorCode MatCompositeSetScalings(Mat,const PetscScalar*);

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

293: PETSC_EXTERN PetscErrorCode MatCreateTranspose(Mat,Mat*);
294: PETSC_EXTERN PetscErrorCode MatTransposeGetMat(Mat,Mat*);
295: PETSC_EXTERN PetscErrorCode MatCreateHermitianTranspose(Mat,Mat*);
296: PETSC_EXTERN PetscErrorCode MatHermitianTransposeGetMat(Mat,Mat*);
297: PETSC_EXTERN PetscErrorCode MatCreateSubMatrixVirtual(Mat,IS,IS,Mat*);
298: PETSC_EXTERN PetscErrorCode MatSubMatrixVirtualUpdate(Mat,Mat,IS,IS);
299: PETSC_EXTERN PetscErrorCode MatCreateLocalRef(Mat,IS,IS,Mat*);
300: PETSC_EXTERN PetscErrorCode MatCreateConstantDiagonal(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar,Mat*);

302: #if defined(PETSC_HAVE_HYPRE)
303: PETSC_EXTERN PetscErrorCode MatHYPRESetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
304: #endif

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

308: PETSC_EXTERN PetscErrorCode MatResetPreallocation(Mat);
309: PETSC_EXTERN PetscErrorCode MatSetUp(Mat);
310: PETSC_EXTERN PetscErrorCode MatDestroy(Mat*);
311: PETSC_EXTERN PetscErrorCode MatGetNonzeroState(Mat,PetscObjectState*);

313: PETSC_EXTERN PetscErrorCode MatConjugate(Mat);
314: PETSC_EXTERN PetscErrorCode MatRealPart(Mat);
315: PETSC_EXTERN PetscErrorCode MatImaginaryPart(Mat);
316: PETSC_EXTERN PetscErrorCode MatGetDiagonalBlock(Mat,Mat*);
317: PETSC_EXTERN PetscErrorCode MatGetTrace(Mat,PetscScalar*);
318: PETSC_EXTERN PetscErrorCode MatInvertBlockDiagonal(Mat,const PetscScalar **);
319: PETSC_EXTERN PetscErrorCode MatInvertVariableBlockDiagonal(Mat,PetscInt,const PetscInt*,PetscScalar*);
320: PETSC_EXTERN PetscErrorCode MatInvertBlockDiagonalMat(Mat,Mat);

322: /* ------------------------------------------------------------*/
323: PETSC_EXTERN PetscErrorCode MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
324: PETSC_EXTERN PetscErrorCode MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
325: PETSC_EXTERN PetscErrorCode MatSetValuesRow(Mat,PetscInt,const PetscScalar[]);
326: PETSC_EXTERN PetscErrorCode MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]);
327: PETSC_EXTERN PetscErrorCode MatSetValuesBatch(Mat,PetscInt,PetscInt,PetscInt[],const PetscScalar[]);
328: PETSC_EXTERN PetscErrorCode MatSetRandom(Mat,PetscRandom);

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

334:    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).
335:    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.

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

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

341:    Level: beginner

343: .seealso:  MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockedStencil(), DMDAVecGetArray(), DMDAVecGetArrayF90()
344: S*/
345: typedef struct {
346:   PetscInt k,j,i,c;
347: } MatStencil;

349: PETSC_EXTERN PetscErrorCode MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
350: PETSC_EXTERN PetscErrorCode MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
351: PETSC_EXTERN PetscErrorCode MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt);

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

357:     Level: beginner

359: .seealso: MatAssemblyBegin(), MatAssemblyEnd()
360: E*/
361: typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType;
362: PETSC_EXTERN PetscErrorCode MatAssemblyBegin(Mat,MatAssemblyType);
363: PETSC_EXTERN PetscErrorCode MatAssemblyEnd(Mat,MatAssemblyType);
364: PETSC_EXTERN PetscErrorCode MatAssembled(Mat,PetscBool *);



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

371:     Level: beginner

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

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

379: .seealso: MatSetOption()
380: E*/
381: typedef enum {MAT_OPTION_MIN = -3,
382:               MAT_UNUSED_NONZERO_LOCATION_ERR = -2,
383:               MAT_ROW_ORIENTED = -1,
384:               MAT_SYMMETRIC = 1,
385:               MAT_STRUCTURALLY_SYMMETRIC = 2,
386:               MAT_NEW_DIAGONALS = 3,
387:               MAT_IGNORE_OFF_PROC_ENTRIES = 4,
388:               MAT_USE_HASH_TABLE = 5,
389:               MAT_KEEP_NONZERO_PATTERN = 6,
390:               MAT_IGNORE_ZERO_ENTRIES = 7,
391:               MAT_USE_INODES = 8,
392:               MAT_HERMITIAN = 9,
393:               MAT_SYMMETRY_ETERNAL = 10,
394:               MAT_NEW_NONZERO_LOCATION_ERR = 11,
395:               MAT_IGNORE_LOWER_TRIANGULAR = 12,
396:               MAT_ERROR_LOWER_TRIANGULAR = 13,
397:               MAT_GETROW_UPPERTRIANGULAR = 14,
398:               MAT_SPD = 15,
399:               MAT_NO_OFF_PROC_ZERO_ROWS = 16,
400:               MAT_NO_OFF_PROC_ENTRIES = 17,
401:               MAT_NEW_NONZERO_LOCATIONS = 18,
402:               MAT_NEW_NONZERO_ALLOCATION_ERR = 19,
403:               MAT_SUBSET_OFF_PROC_ENTRIES = 20,
404:               MAT_SUBMAT_SINGLEIS = 21,
405:               MAT_STRUCTURE_ONLY = 22,
406:               MAT_SORTED_FULL = 23,
407:               MAT_OPTION_MAX = 24} MatOption;

409: PETSC_EXTERN const char *const *MatOptions;
410: PETSC_EXTERN PetscErrorCode MatSetOption(Mat,MatOption,PetscBool);
411: PETSC_EXTERN PetscErrorCode MatGetOption(Mat,MatOption,PetscBool*);
412: PETSC_EXTERN PetscErrorCode MatGetType(Mat,MatType*);

414: PETSC_EXTERN PetscErrorCode MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
415: PETSC_EXTERN PetscErrorCode MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
416: PETSC_EXTERN PetscErrorCode MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
417: PETSC_EXTERN PetscErrorCode MatGetRowUpperTriangular(Mat);
418: PETSC_EXTERN PetscErrorCode MatRestoreRowUpperTriangular(Mat);
419: PETSC_EXTERN PetscErrorCode MatGetColumnVector(Mat,Vec,PetscInt);
420: PETSC_EXTERN PetscErrorCode MatSeqAIJGetArray(Mat,PetscScalar *[]);
421: PETSC_EXTERN PetscErrorCode MatSeqAIJRestoreArray(Mat,PetscScalar *[]);
422: PETSC_EXTERN PetscErrorCode MatSeqAIJGetMaxRowNonzeros(Mat,PetscInt*);
423: PETSC_EXTERN PetscErrorCode MatSeqAIJSetValuesLocalFast(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
424: PETSC_EXTERN PetscErrorCode MatSeqAIJSetType(Mat,MatType);
425: PETSC_EXTERN PetscErrorCode MatSeqAIJRegister(const char[],PetscErrorCode (*)(Mat,MatType,MatReuse,Mat *));
426: PETSC_EXTERN PetscFunctionList MatSeqAIJList;
427: PETSC_EXTERN PetscErrorCode MatSeqSBAIJGetArray(Mat,PetscScalar *[]);
428: PETSC_EXTERN PetscErrorCode MatSeqSBAIJRestoreArray(Mat,PetscScalar *[]);
429: PETSC_EXTERN PetscErrorCode MatDenseGetArray(Mat,PetscScalar *[]);
430: PETSC_EXTERN PetscErrorCode MatDenseRestoreArray(Mat,PetscScalar *[]);
431: PETSC_EXTERN PetscErrorCode MatDensePlaceArray(Mat,const PetscScalar[]);
432: PETSC_EXTERN PetscErrorCode MatDenseResetArray(Mat);
433: PETSC_EXTERN PetscErrorCode MatDenseGetArrayRead(Mat,const PetscScalar *[]);
434: PETSC_EXTERN PetscErrorCode MatDenseRestoreArrayRead(Mat,const PetscScalar *[]);
435: PETSC_EXTERN PetscErrorCode MatGetBlockSize(Mat,PetscInt *);
436: PETSC_EXTERN PetscErrorCode MatSetBlockSize(Mat,PetscInt);
437: PETSC_EXTERN PetscErrorCode MatGetBlockSizes(Mat,PetscInt *,PetscInt *);
438: PETSC_EXTERN PetscErrorCode MatSetBlockSizes(Mat,PetscInt,PetscInt);
439: PETSC_EXTERN PetscErrorCode MatSetBlockSizesFromMats(Mat,Mat,Mat);
440: PETSC_EXTERN PetscErrorCode MatSetVariableBlockSizes(Mat,PetscInt,PetscInt*);
441: PETSC_EXTERN PetscErrorCode MatGetVariableBlockSizes(Mat,PetscInt*,const PetscInt**);

443: PETSC_EXTERN PetscErrorCode MatDenseGetColumn(Mat,PetscInt,PetscScalar *[]);
444: PETSC_EXTERN PetscErrorCode MatDenseRestoreColumn(Mat,PetscScalar *[]);

446: PETSC_EXTERN PetscErrorCode MatMult(Mat,Vec,Vec);
447: PETSC_EXTERN PetscErrorCode MatMultDiagonalBlock(Mat,Vec,Vec);
448: PETSC_EXTERN PetscErrorCode MatMultAdd(Mat,Vec,Vec,Vec);
449: PETSC_EXTERN PetscErrorCode MatMultTranspose(Mat,Vec,Vec);
450: PETSC_EXTERN PetscErrorCode MatMultHermitianTranspose(Mat,Vec,Vec);
451: PETSC_EXTERN PetscErrorCode MatIsTranspose(Mat,Mat,PetscReal,PetscBool *);
452: PETSC_EXTERN PetscErrorCode MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *);
453: PETSC_EXTERN PetscErrorCode MatMultTransposeAdd(Mat,Vec,Vec,Vec);
454: PETSC_EXTERN PetscErrorCode MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec);
455: PETSC_EXTERN PetscErrorCode MatMultConstrained(Mat,Vec,Vec);
456: PETSC_EXTERN PetscErrorCode MatMultTransposeConstrained(Mat,Vec,Vec);
457: PETSC_EXTERN PetscErrorCode MatMatSolve(Mat,Mat,Mat);
458: PETSC_EXTERN PetscErrorCode MatMatSolveTranspose(Mat,Mat,Mat);
459: PETSC_EXTERN PetscErrorCode MatMatTransposeSolve(Mat,Mat,Mat);
460: PETSC_EXTERN PetscErrorCode MatResidual(Mat,Vec,Vec,Vec);

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

466:     Level: beginner

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

470: $   MAT_DO_NOT_COPY_VALUES    - Create a matrix using the same nonzero pattern as the original matrix,
471: $                               with zeros for the numerical values.
472: $   MAT_COPY_VALUES           - Create a matrix with the same nonzero pattern as the original matrix
473: $                               and with the same numerical values.
474: $   MAT_SHARE_NONZERO_PATTERN - Create a matrix that shares the nonzero structure with the previous matrix
475: $                               and does not copy it, using zeros for the numerical values. The parent and
476: $                               child matrices will share their index (i and j) arrays, and you cannot
477: $                               insert new nonzero entries into either matrix.

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

483: .seealso: MatDuplicate()
484: E*/
485: typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption;

487: PETSC_EXTERN PetscErrorCode MatConvert(Mat,MatType,MatReuse,Mat*);
488: PETSC_EXTERN PetscErrorCode MatDuplicate(Mat,MatDuplicateOption,Mat*);


491: PETSC_EXTERN PetscErrorCode MatCopy(Mat,Mat,MatStructure);
492: PETSC_EXTERN PetscErrorCode MatView(Mat,PetscViewer);
493: PETSC_EXTERN PetscErrorCode MatIsSymmetric(Mat,PetscReal,PetscBool *);
494: PETSC_EXTERN PetscErrorCode MatIsStructurallySymmetric(Mat,PetscBool *);
495: PETSC_EXTERN PetscErrorCode MatIsHermitian(Mat,PetscReal,PetscBool *);
496: PETSC_EXTERN PetscErrorCode MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *);
497: PETSC_EXTERN PetscErrorCode MatIsHermitianKnown(Mat,PetscBool *,PetscBool *);
498: PETSC_EXTERN PetscErrorCode MatMissingDiagonal(Mat,PetscBool  *,PetscInt *);
499: PETSC_EXTERN PetscErrorCode MatLoad(Mat, PetscViewer);

501: PETSC_EXTERN PetscErrorCode MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
502: PETSC_EXTERN PetscErrorCode MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool  *);
503: PETSC_EXTERN PetscErrorCode MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
504: PETSC_EXTERN PetscErrorCode MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool  *);

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

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

511:    Level: intermediate

513: .seealso:  MatGetInfo(), MatInfoType
514: S*/
515: typedef struct {
516:   PetscLogDouble block_size;                         /* block size */
517:   PetscLogDouble nz_allocated,nz_used,nz_unneeded;   /* number of nonzeros */
518:   PetscLogDouble memory;                             /* memory allocated */
519:   PetscLogDouble assemblies;                         /* number of matrix assemblies called */
520:   PetscLogDouble mallocs;                            /* number of mallocs during MatSetValues() */
521:   PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
522:   PetscLogDouble factor_mallocs;                     /* number of mallocs during factorization */
523: } MatInfo;

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

529:     Level: beginner

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

533: .seealso: MatGetInfo(), MatInfo
534: E*/
535: typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
536: PETSC_EXTERN PetscErrorCode MatGetInfo(Mat,MatInfoType,MatInfo*);
537: PETSC_EXTERN PetscErrorCode MatGetDiagonal(Mat,Vec);
538: PETSC_EXTERN PetscErrorCode MatGetRowMax(Mat,Vec,PetscInt[]);
539: PETSC_EXTERN PetscErrorCode MatGetRowMin(Mat,Vec,PetscInt[]);
540: PETSC_EXTERN PetscErrorCode MatGetRowMaxAbs(Mat,Vec,PetscInt[]);
541: PETSC_EXTERN PetscErrorCode MatGetRowMinAbs(Mat,Vec,PetscInt[]);
542: PETSC_EXTERN PetscErrorCode MatGetRowSum(Mat,Vec);
543: PETSC_EXTERN PetscErrorCode MatTranspose(Mat,MatReuse,Mat*);
544: PETSC_EXTERN PetscErrorCode MatHermitianTranspose(Mat,MatReuse,Mat*);
545: PETSC_EXTERN PetscErrorCode MatPermute(Mat,IS,IS,Mat*);
546: PETSC_EXTERN PetscErrorCode MatDiagonalScale(Mat,Vec,Vec);
547: PETSC_EXTERN PetscErrorCode MatDiagonalSet(Mat,Vec,InsertMode);

549: PETSC_EXTERN PetscErrorCode MatEqual(Mat,Mat,PetscBool*);
550: PETSC_EXTERN PetscErrorCode MatMultEqual(Mat,Mat,PetscInt,PetscBool*);
551: PETSC_EXTERN PetscErrorCode MatMultAddEqual(Mat,Mat,PetscInt,PetscBool*);
552: PETSC_EXTERN PetscErrorCode MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool*);
553: PETSC_EXTERN PetscErrorCode MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool*);
554: PETSC_EXTERN PetscErrorCode MatMatMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*);
555: PETSC_EXTERN PetscErrorCode MatTransposeMatMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*);
556: PETSC_EXTERN PetscErrorCode MatMatTransposeMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*);
557: PETSC_EXTERN PetscErrorCode MatIsLinear(Mat,PetscInt,PetscBool*);

559: PETSC_EXTERN PetscErrorCode MatNorm(Mat,NormType,PetscReal*);
560: PETSC_EXTERN PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal*);
561: PETSC_EXTERN PetscErrorCode MatZeroEntries(Mat);
562: PETSC_EXTERN PetscErrorCode MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
563: PETSC_EXTERN PetscErrorCode MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec);
564: PETSC_EXTERN PetscErrorCode MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec);
565: PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsStencil(Mat,PetscInt,const MatStencil[],PetscScalar,Vec,Vec);
566: PETSC_EXTERN PetscErrorCode MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
567: PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec);

569: PETSC_EXTERN PetscErrorCode MatGetSize(Mat,PetscInt*,PetscInt*);
570: PETSC_EXTERN PetscErrorCode MatGetLocalSize(Mat,PetscInt*,PetscInt*);
571: PETSC_EXTERN PetscErrorCode MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
572: PETSC_EXTERN PetscErrorCode MatGetOwnershipRanges(Mat,const PetscInt**);
573: PETSC_EXTERN PetscErrorCode MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*);
574: PETSC_EXTERN PetscErrorCode MatGetOwnershipRangesColumn(Mat,const PetscInt**);
575: PETSC_EXTERN PetscErrorCode MatGetOwnershipIS(Mat,IS*,IS*);

577: PETSC_EXTERN PetscErrorCode MatCreateSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
578: 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);}
579: PETSC_EXTERN PetscErrorCode MatCreateSubMatricesMPI(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
580: 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);}
581: PETSC_EXTERN PetscErrorCode MatDestroyMatrices(PetscInt,Mat *[]);
582: PETSC_EXTERN PetscErrorCode MatDestroySubMatrices(PetscInt,Mat *[]);
583: PETSC_EXTERN PetscErrorCode MatCreateSubMatrix(Mat,IS,IS,MatReuse,Mat *);
584: 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);}
585: PETSC_EXTERN PetscErrorCode MatGetLocalSubMatrix(Mat,IS,IS,Mat*);
586: PETSC_EXTERN PetscErrorCode MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*);
587: PETSC_EXTERN PetscErrorCode MatGetSeqNonzeroStructure(Mat,Mat*);
588: PETSC_EXTERN PetscErrorCode MatDestroySeqNonzeroStructure(Mat*);

590: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJ(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*);
591: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*);
592: PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJNumeric(Mat,Mat);
593: PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMat(Mat,MatReuse,Mat*);
594: PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*);
595: PETSC_EXTERN PetscErrorCode MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,Mat*);
596: PETSC_EXTERN PetscErrorCode MatGetGhosts(Mat, PetscInt *,const PetscInt *[]);

598: PETSC_EXTERN PetscErrorCode MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt);
599: PETSC_EXTERN PetscErrorCode MatIncreaseOverlapSplit(Mat mat,PetscInt n,IS is[],PetscInt ov);
600: PETSC_EXTERN PetscErrorCode MatMPIAIJSetUseScalableIncreaseOverlap(Mat,PetscBool);

602: PETSC_EXTERN PetscErrorCode MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
603: PETSC_EXTERN PetscErrorCode MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
604: PETSC_EXTERN PetscErrorCode MatMatMultNumeric(Mat,Mat,Mat);

606: PETSC_EXTERN PetscErrorCode MatMatMatMult(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
607: PETSC_EXTERN PetscErrorCode MatGalerkin(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);

609: PETSC_EXTERN PetscErrorCode MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*);
610: PETSC_EXTERN PetscErrorCode MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*);
611: PETSC_EXTERN PetscErrorCode MatPtAPNumeric(Mat,Mat,Mat);
612: PETSC_EXTERN PetscErrorCode MatRARt(Mat,Mat,MatReuse,PetscReal,Mat*);
613: PETSC_EXTERN PetscErrorCode MatRARtSymbolic(Mat,Mat,PetscReal,Mat*);
614: PETSC_EXTERN PetscErrorCode MatRARtNumeric(Mat,Mat,Mat);

616: PETSC_EXTERN PetscErrorCode MatTransposeMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
617: PETSC_EXTERN PetscErrorCode MatTransposetMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
618: PETSC_EXTERN PetscErrorCode MatTransposetMatMultNumeric(Mat,Mat,Mat);
619: PETSC_EXTERN PetscErrorCode MatMatTransposeMult(Mat,Mat,MatReuse,PetscReal,Mat*);
620: PETSC_EXTERN PetscErrorCode MatMatTransposeMultSymbolic(Mat,Mat,PetscReal,Mat*);
621: PETSC_EXTERN PetscErrorCode MatMatTransposeMultNumeric(Mat,Mat,Mat);

623: PETSC_EXTERN PetscErrorCode MatAXPY(Mat,PetscScalar,Mat,MatStructure);
624: PETSC_EXTERN PetscErrorCode MatAYPX(Mat,PetscScalar,Mat,MatStructure);

626: PETSC_EXTERN PetscErrorCode MatScale(Mat,PetscScalar);
627: PETSC_EXTERN PetscErrorCode MatShift(Mat,PetscScalar);

629: PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
630: PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*);
631: PETSC_EXTERN PetscErrorCode MatGetLayouts(Mat,PetscLayout*,PetscLayout*);
632: PETSC_EXTERN PetscErrorCode MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
633: PETSC_EXTERN PetscErrorCode MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
634: PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
635: PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
636: PETSC_EXTERN PetscErrorCode MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
637: PETSC_EXTERN PetscErrorCode MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);

639: PETSC_EXTERN PetscErrorCode MatStashSetInitialSize(Mat,PetscInt,PetscInt);
640: PETSC_EXTERN PetscErrorCode MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);

642: PETSC_EXTERN PetscErrorCode MatInterpolate(Mat,Vec,Vec);
643: PETSC_EXTERN PetscErrorCode MatInterpolateAdd(Mat,Vec,Vec,Vec);
644: PETSC_EXTERN PetscErrorCode MatRestrict(Mat,Vec,Vec);
645: PETSC_EXTERN PetscErrorCode MatCreateVecs(Mat,Vec*,Vec*);
646: 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);}
647: PETSC_EXTERN PetscErrorCode MatCreateRedundantMatrix(Mat,PetscInt,MPI_Comm,MatReuse,Mat*);
648: PETSC_EXTERN PetscErrorCode MatGetMultiProcBlock(Mat,MPI_Comm,MatReuse,Mat*);
649: PETSC_EXTERN PetscErrorCode MatFindZeroDiagonals(Mat,IS*);
650: PETSC_EXTERN PetscErrorCode MatFindOffBlockDiagonalEntries(Mat,IS*);
651: PETSC_EXTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);

653: /*MC
654:    MatSetValue - Set a single entry into a matrix.

656:    Not collective

658:    Synopsis:
659:      #include <petscmat.h>
660:      PetscErrorCode MatSetValue(Mat m,PetscInt row,PetscInt col,PetscScalar value,InsertMode mode)

662:    Input Parameters:
663: +  m - the matrix
664: .  row - the row location of the entry
665: .  col - the column location of the entry
666: .  value - the value to insert
667: -  mode - either INSERT_VALUES or ADD_VALUES

669:    Notes:
670:    For efficiency one should use MatSetValues() and set several or many
671:    values simultaneously if possible.

673:    Level: beginner

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

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

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

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

687:    Synopsis:
688:    #include <petscmat.h>
689:    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)

691:    Collective

693:    Input Parameters:
694: +  comm - the communicator that will share the eventually allocated matrix
695: .  nrows - the number of LOCAL rows in the matrix
696: -  ncols - the number of LOCAL columns in the matrix

698:    Output Parameters:
699: +  dnz - the array that will be passed to the matrix preallocation routines
700: -  onz - the other array passed to the matrix preallocation routines

702:    Level: intermediate

704:    Notes:
705:     See Users-Manual: Chapter 14 Hints for Performance Tuning for more details.

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

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

711: .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
712:           MatPreallocateSymmetricSetLocalBlock()
713: M*/
714: #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
715: do { \
716:   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ncols = (ncols),__rstart,__start,__end; \
717:   _4_PetscCalloc2((size_t)__nrows,&dnz,(size_t)__nrows,&onz);CHKERRQ(_4_ierr); \
718:   __start = 0; __end = __start;                                         \
719:   _4_MPI_Scan(&__ncols,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ncols;\
720:   _4_MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows; \
721:   do { } while(0)

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

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

731:    Not Collective

733:    Input Parameters:
734: +  map - the row mapping from local numbering to global numbering
735: .  nrows - the number of rows indicated
736: .  rows - the indices of the rows
737: .  cmap - the column mapping from local to global numbering
738: .  ncols - the number of columns in the matrix
739: .  cols - the columns indicated
740: .  dnz - the array that will be passed to the matrix preallocation routines
741: -  onz - the other array passed to the matrix preallocation routines

743:    Level: intermediate

745:    Notes:
746:     See Users-Manual: Chapter 14 Hints for Performance Tuning for more details.

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

750: .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock()
751:           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocalRemoveDups()
752: M*/
753: #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
754: do {\
755:   PetscInt __l;\
756:   _4_ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
757:   _4_ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
758:   for (__l=0;__l<nrows;__l++) {\
759:     _4_MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
760:   }\
761:  } while (0)

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

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

771:    Not Collective

773:    Input Parameters:
774: +  map - the row mapping from local numbering to global numbering
775: .  nrows - the number of rows indicated
776: .  rows - the indices of the rows (these values are mapped to the global values)
777: .  cmap - the column mapping from local to global numbering
778: .  ncols - the number of columns in the matrix   (this value will be changed if duplicate columns are found)
779: .  cols - the columns indicated (these values are mapped to the global values, they are then sorted and duplicates removed)
780: .  dnz - the array that will be passed to the matrix preallocation routines
781: -  onz - the other array passed to the matrix preallocation routines

783:    Level: intermediate

785:    Notes:
786:     See Users-Manual: Chapter 14 Hints for Performance Tuning for more details.

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

790: .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock()
791:           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal()
792: M*/
793: #define MatPreallocateSetLocalRemoveDups(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
794: do {\
795:   PetscInt __l;\
796:   _4_ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
797:   _4_ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
798:   _4_PetscSortRemoveDupsInt(&ncols,cols);CHKERRQ(_4_ierr);\
799:   for (__l=0;__l<nrows;__l++) {\
800:     _4_MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
801:   }\
802:  } while (0)

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

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

812:    Not Collective

814:    Input Parameters:
815: +  map - the row mapping from local numbering to global numbering
816: .  nrows - the number of rows indicated
817: .  rows - the indices of the rows
818: .  cmap - the column mapping from local to global numbering
819: .  ncols - the number of columns in the matrix
820: .  cols - the columns indicated
821: .  dnz - the array that will be passed to the matrix preallocation routines
822: -  onz - the other array passed to the matrix preallocation routines

824:    Level: intermediate

826:    Notes:
827:     See Users-Manual: Chapter 14 Hints for Performance Tuning for more details.

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

831: .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock()
832:           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock()
833: M*/
834: #define MatPreallocateSetLocalBlock(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
835: do {\
836:   PetscInt __l;\
837:   _4_ISLocalToGlobalMappingApplyBlock(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
838:   _4_ISLocalToGlobalMappingApplyBlock(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
839:   for (__l=0;__l<nrows;__l++) {\
840:     _4_MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
841:   }\
842: } while (0)

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

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

852:    Not Collective

854:    Input Parameters:
855: +  map - the mapping between local numbering and global numbering
856: .  nrows - the number of rows indicated
857: .  rows - the indices of the rows
858: .  ncols - the number of columns in the matrix
859: .  cols - the columns indicated
860: .  dnz - the array that will be passed to the matrix preallocation routines
861: -  onz - the other array passed to the matrix preallocation routines

863:    Level: intermediate

865:    Notes:
866:     See Users-Manual: Chapter 14 Hints for Performance Tuning for more details.

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

870: .seealso: MatPreallocateFinalize(), MatPreallocateSet()
871:           MatPreallocateInitialize(),  MatPreallocateSetLocal()
872: M*/
873: #define MatPreallocateSymmetricSetLocalBlock(map,nrows,rows,ncols,cols,dnz,onz) 0;\
874: do {\
875:   PetscInt __l;\
876:   _4_ISLocalToGlobalMappingApplyBlock(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
877:   _4_ISLocalToGlobalMappingApplyBlock(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
878:   for (__l=0;__l<nrows;__l++) {\
879:     _4_MatPreallocateSymmetricSetBlock((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
880:   }\
881: } while (0)

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

887:    Synopsis:
888:    #include <petscmat.h>
889:    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)

891:    Not Collective

893:    Input Parameters:
894: +  row - the row
895: .  ncols - the number of columns in the matrix
896: -  cols - the columns indicated

898:    Output Parameters:
899: +  dnz - the array that will be passed to the matrix preallocation routines
900: -  onz - the other array passed to the matrix preallocation routines

902:    Level: intermediate

904:    Notes:
905:     See Users-Manual: Chapter 14 Hints for Performance Tuning for more details.

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

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

911: .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock()
912:           MatPreallocateInitialize(), MatPreallocateSetLocal()
913: M*/
914: #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
915: do { PetscInt __i; \
916:   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);\
917:   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);\
918:   for (__i=0; __i<nc; __i++) {\
919:     if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \
920:     else if (dnz[row - __rstart] < __ncols) dnz[row - __rstart]++;\
921:   }\
922: } while (0)

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

928:    Synopsis:
929:    #include <petscmat.h>
930:    PetscErrorCode MatPreallocateSymmetricSetBlock(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)

932:    Not Collective

934:    Input Parameters:
935: +  nrows - the number of rows indicated
936: .  rows - the indices of the rows
937: .  ncols - the number of columns in the matrix
938: .  cols - the columns indicated
939: .  dnz - the array that will be passed to the matrix preallocation routines
940: -  onz - the other array passed to the matrix preallocation routines

942:    Level: intermediate

944:    Notes:
945:     See Users-Manual: Chapter 14 Hints for Performance Tuning for more details.

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

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

951: .seealso: MatPreallocateFinalize(), MatPreallocateSet(),  MatPreallocateInitialize(),
952:           MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal()
953: M*/
954: #define MatPreallocateSymmetricSetBlock(row,nc,cols,dnz,onz) 0;\
955: do { PetscInt __i; \
956:   for (__i=0; __i<nc; __i++) {\
957:     if (cols[__i] >= __end) onz[row - __rstart]++; \
958:     else if (cols[__i] >= row && dnz[row - __rstart] < __ncols) dnz[row - __rstart]++;\
959:   }\
960: } while (0)

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

965:    Synopsis:
966:    #include <petscmat.h>
967:    PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz)

969:    Not Collective

971:    Input Parameters:
972: +  A - matrix
973: .  row - row where values exist (must be local to this process)
974: .  ncols - number of columns
975: .  cols - columns with nonzeros
976: .  dnz - the array that will be passed to the matrix preallocation routines
977: -  onz - the other array passed to the matrix preallocation routines

979:    Level: intermediate

981:    Notes:
982:     See Users-Manual: Chapter 14 Hints for Performance Tuning for more details.

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

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

988: .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
989:           MatPreallocateSymmetricSetLocalBlock()
990: M*/
991: #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)


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

998:    Synopsis:
999:    #include <petscmat.h>
1000:    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)

1002:    Collective

1004:    Input Parameters:
1005: +  dnz - the array that was be passed to the matrix preallocation routines
1006: -  onz - the other array passed to the matrix preallocation routines

1008:    Level: intermediate

1010:    Notes:
1011:     See Users-Manual: Chapter 14 Hints for Performance Tuning for more details.

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

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

1017: .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
1018:           MatPreallocateSymmetricSetLocalBlock()
1019: M*/
1020: #define MatPreallocateFinalize(dnz,onz) 0;_4_PetscFree2(dnz,onz);CHKERRQ(_4_ierr);} while(0)

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

1025: PETSC_EXTERN PetscErrorCode MatInodeAdjustForInodes(Mat,IS*,IS*);
1026: PETSC_EXTERN PetscErrorCode MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);

1028: PETSC_EXTERN PetscErrorCode MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
1029: PETSC_EXTERN PetscErrorCode MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
1030: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
1031: PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
1032: PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
1033: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJFromTriple(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*,PetscInt,PetscBool);

1035: #define MAT_SKIP_ALLOCATION -4

1037: PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
1038: PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
1039: PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);

1041: PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1042: PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1043: PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1044: PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []);
1045: PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
1046: PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
1047: PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
1048: PETSC_EXTERN PetscErrorCode MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
1049: PETSC_EXTERN PetscErrorCode MatMPIAdjToSeq(Mat,Mat*);
1050: PETSC_EXTERN PetscErrorCode MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
1051: PETSC_EXTERN PetscErrorCode MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
1052: PETSC_EXTERN PetscErrorCode MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,const PetscInt*[]);
1053: PETSC_EXTERN PetscErrorCode MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,const PetscInt*[]);
1054: PETSC_EXTERN PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat,Mat*);

1056: PETSC_EXTERN PetscErrorCode MatSeqDenseSetLDA(Mat,PetscInt);
1057: PETSC_EXTERN PetscErrorCode MatDenseGetLDA(Mat,PetscInt*);
1058: PETSC_EXTERN PetscErrorCode MatDenseGetLocalMatrix(Mat,Mat*);

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

1062: PETSC_EXTERN PetscErrorCode MatStoreValues(Mat);
1063: PETSC_EXTERN PetscErrorCode MatRetrieveValues(Mat);

1065: PETSC_EXTERN PetscErrorCode MatDAADSetCtx(Mat,void*);

1067: PETSC_EXTERN PetscErrorCode MatFindNonzeroRows(Mat,IS*);
1068: PETSC_EXTERN PetscErrorCode MatFindZeroRows(Mat,IS*);
1069: /*
1070:   These routines are not usually accessed directly, rather solving is
1071:   done through the KSP and PC interfaces.
1072: */

1074: /*J
1075:     MatOrderingType - String with the name of a PETSc matrix ordering

1077:    Level: beginner

1079: .seealso: MatGetOrdering()
1080: J*/
1081: typedef const char* MatOrderingType;
1082: #define MATORDERINGNATURAL     "natural"
1083: #define MATORDERINGND          "nd"
1084: #define MATORDERING1WD         "1wd"
1085: #define MATORDERINGRCM         "rcm"
1086: #define MATORDERINGQMD         "qmd"
1087: #define MATORDERINGROWLENGTH   "rowlength"
1088: #define MATORDERINGWBM         "wbm"
1089: #define MATORDERINGSPECTRAL    "spectral"
1090: #define MATORDERINGAMD         "amd"            /* only works if UMFPACK is installed with PETSc */

1092: PETSC_EXTERN PetscErrorCode MatGetOrdering(Mat,MatOrderingType,IS*,IS*);
1093: PETSC_EXTERN PetscErrorCode MatGetOrderingList(PetscFunctionList*);
1094: PETSC_EXTERN PetscErrorCode MatOrderingRegister(const char[],PetscErrorCode(*)(Mat,MatOrderingType,IS*,IS*));
1095: PETSC_EXTERN PetscFunctionList MatOrderingList;

1097: PETSC_EXTERN PetscErrorCode MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
1098: PETSC_EXTERN PetscErrorCode MatCreateLaplacian(Mat,PetscReal,PetscBool,Mat*);

1100: /*S
1101:     MatFactorShiftType - Numeric Shift.

1103:    Level: beginner

1105: S*/
1106: typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType;
1107: PETSC_EXTERN const char *const MatFactorShiftTypes[];
1108: PETSC_EXTERN const char *const MatFactorShiftTypesDetail[];

1110: /*S
1111:     MatFactorError - indicates what type of error in matrix factor

1113:     Level: beginner

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

1118: .seealso: MatGetFactor()
1119: S*/
1120: typedef enum {MAT_FACTOR_NOERROR,MAT_FACTOR_STRUCT_ZEROPIVOT,MAT_FACTOR_NUMERIC_ZEROPIVOT,MAT_FACTOR_OUTMEMORY,MAT_FACTOR_OTHER} MatFactorError;

1122: PETSC_EXTERN PetscErrorCode MatFactorGetError(Mat,MatFactorError*);
1123: PETSC_EXTERN PetscErrorCode MatFactorClearError(Mat);
1124: PETSC_EXTERN PetscErrorCode MatFactorGetErrorZeroPivot(Mat,PetscReal*,PetscInt*);

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

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

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

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

1137:    Level: developer

1139: .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
1140:           MatFactorInfoInitialize()

1142: S*/
1143: typedef struct {
1144:   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
1145:   PetscReal     usedt;
1146:   PetscReal     dt;             /* drop tolerance */
1147:   PetscReal     dtcol;          /* tolerance for pivoting */
1148:   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
1149:   PetscReal     fill;           /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */
1150:   PetscReal     levels;         /* ICC/ILU(levels) */
1151:   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
1152:                                    factorization may be faster if do not pivot */
1153:   PetscReal     zeropivot;      /* pivot is called zero if less than this */
1154:   PetscReal     shifttype;      /* type of shift added to matrix factor to prevent zero pivots */
1155:   PetscReal     shiftamount;     /* how large the shift is */
1156: } MatFactorInfo;

1158: PETSC_EXTERN PetscErrorCode MatFactorInfoInitialize(MatFactorInfo*);
1159: PETSC_EXTERN PetscErrorCode MatCholeskyFactor(Mat,IS,const MatFactorInfo*);
1160: PETSC_EXTERN PetscErrorCode MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
1161: PETSC_EXTERN PetscErrorCode MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*);
1162: PETSC_EXTERN PetscErrorCode MatLUFactor(Mat,IS,IS,const MatFactorInfo*);
1163: PETSC_EXTERN PetscErrorCode MatILUFactor(Mat,IS,IS,const MatFactorInfo*);
1164: PETSC_EXTERN PetscErrorCode MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
1165: PETSC_EXTERN PetscErrorCode MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
1166: PETSC_EXTERN PetscErrorCode MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
1167: PETSC_EXTERN PetscErrorCode MatICCFactor(Mat,IS,const MatFactorInfo*);
1168: PETSC_EXTERN PetscErrorCode MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*);
1169: PETSC_EXTERN PetscErrorCode MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
1170: PETSC_EXTERN PetscErrorCode MatSolve(Mat,Vec,Vec);
1171: PETSC_EXTERN PetscErrorCode MatForwardSolve(Mat,Vec,Vec);
1172: PETSC_EXTERN PetscErrorCode MatBackwardSolve(Mat,Vec,Vec);
1173: PETSC_EXTERN PetscErrorCode MatSolveAdd(Mat,Vec,Vec,Vec);
1174: PETSC_EXTERN PetscErrorCode MatSolveTranspose(Mat,Vec,Vec);
1175: PETSC_EXTERN PetscErrorCode MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
1176: PETSC_EXTERN PetscErrorCode MatSolves(Mat,Vecs,Vecs);
1177: PETSC_EXTERN PetscErrorCode MatSetUnfactored(Mat);

1179: typedef enum {MAT_FACTOR_SCHUR_UNFACTORED, MAT_FACTOR_SCHUR_FACTORED, MAT_FACTOR_SCHUR_INVERTED} MatFactorSchurStatus;
1180: PETSC_EXTERN PetscErrorCode MatFactorSetSchurIS(Mat,IS);
1181: PETSC_EXTERN PetscErrorCode MatFactorGetSchurComplement(Mat,Mat*,MatFactorSchurStatus*);
1182: PETSC_EXTERN PetscErrorCode MatFactorRestoreSchurComplement(Mat,Mat*,MatFactorSchurStatus);
1183: PETSC_EXTERN PetscErrorCode MatFactorInvertSchurComplement(Mat);
1184: PETSC_EXTERN PetscErrorCode MatFactorCreateSchurComplement(Mat,Mat*,MatFactorSchurStatus*);
1185: PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplement(Mat,Vec,Vec);
1186: PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplementTranspose(Mat,Vec,Vec);
1187: PETSC_EXTERN PetscErrorCode MatFactorFactorizeSchurComplement(Mat);

1189: /*E
1190:     MatSORType - What type of (S)SOR to perform

1192:     Level: beginner

1194:    May be bitwise ORd together

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

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

1200: .seealso: MatSOR()
1201: E*/
1202: typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1203:               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1204:               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
1205:               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
1206: PETSC_EXTERN PetscErrorCode MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);

1208: /*
1209:     These routines are for efficiently computing Jacobians via finite differences.
1210: */

1212: /*S
1213:      MatColoring - Object for managing the coloring of matrices.

1215:    Level: beginner

1217: .seealso:  MatFDColoringCreate() ISColoring MatFDColoring
1218: S*/
1219: typedef struct _p_MatColoring* MatColoring;
1220: /*J
1221:     MatColoringType - String with the name of a PETSc matrix coloring

1223:    Level: beginner

1225: .seealso: MatColoringSetType(), MatColoring
1226: J*/

1228: typedef const  char*           MatColoringType;
1229: #define MATCOLORINGJP      "jp"
1230: #define MATCOLORINGPOWER   "power"
1231: #define MATCOLORINGNATURAL "natural"
1232: #define MATCOLORINGSL      "sl"
1233: #define MATCOLORINGLF      "lf"
1234: #define MATCOLORINGID      "id"
1235: #define MATCOLORINGGREEDY  "greedy"

1237: /*E
1238:    MatColoringWeightType - Type of weight scheme

1240:     Not Collective

1242: +   MAT_COLORING_RANDOM  - Random weights
1243: .   MAT_COLORING_LEXICAL - Lexical weighting based upon global numbering.
1244: -   MAT_COLORING_LF      - Last-first weighting.

1246:     Level: intermediate

1248:    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1249: E*/
1250: typedef enum {MAT_COLORING_WEIGHT_RANDOM,MAT_COLORING_WEIGHT_LEXICAL,MAT_COLORING_WEIGHT_LF,MAT_COLORING_WEIGHT_SL} MatColoringWeightType;

1252: PETSC_EXTERN PetscErrorCode MatColoringCreate(Mat,MatColoring*);
1253: PETSC_EXTERN PetscErrorCode MatColoringGetDegrees(Mat,PetscInt,PetscInt*);
1254: PETSC_EXTERN PetscErrorCode MatColoringDestroy(MatColoring*);
1255: PETSC_EXTERN PetscErrorCode MatColoringView(MatColoring,PetscViewer);
1256: PETSC_EXTERN PetscErrorCode MatColoringSetType(MatColoring,MatColoringType);
1257: PETSC_EXTERN PetscErrorCode MatColoringSetFromOptions(MatColoring);
1258: PETSC_EXTERN PetscErrorCode MatColoringSetDistance(MatColoring,PetscInt);
1259: PETSC_EXTERN PetscErrorCode MatColoringGetDistance(MatColoring,PetscInt*);
1260: PETSC_EXTERN PetscErrorCode MatColoringSetMaxColors(MatColoring,PetscInt);
1261: PETSC_EXTERN PetscErrorCode MatColoringGetMaxColors(MatColoring,PetscInt*);
1262: PETSC_EXTERN PetscErrorCode MatColoringApply(MatColoring,ISColoring*);
1263: PETSC_EXTERN PetscErrorCode MatColoringRegister(const char[],PetscErrorCode(*)(MatColoring));
1264: PETSC_EXTERN PetscErrorCode MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
1265: PETSC_EXTERN PetscErrorCode MatColoringSetWeightType(MatColoring,MatColoringWeightType);
1266: PETSC_EXTERN PetscErrorCode MatColoringSetWeights(MatColoring,PetscReal*,PetscInt*);
1267: PETSC_EXTERN PetscErrorCode MatColoringCreateWeights(MatColoring,PetscReal **,PetscInt **lperm);
1268: PETSC_EXTERN PetscErrorCode MatColoringTest(MatColoring,ISColoring);
1269: PETSC_DEPRECATED_FUNCTION("Use MatColoringTest() (since version 3.10)") PETSC_STATIC_INLINE PetscErrorCode MatColoringTestValid(MatColoring matcoloring,ISColoring iscoloring) {return MatColoringTest(matcoloring,iscoloring);}
1270: PETSC_EXTERN PetscErrorCode MatISColoringTest(Mat,ISColoring);

1272: /*S
1273:      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1274:         and coloring

1276:    Level: beginner

1278: .seealso:  MatFDColoringCreate()
1279: S*/
1280: typedef struct _p_MatFDColoring* MatFDColoring;

1282: PETSC_EXTERN PetscErrorCode MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
1283: PETSC_EXTERN PetscErrorCode MatFDColoringDestroy(MatFDColoring*);
1284: PETSC_EXTERN PetscErrorCode MatFDColoringView(MatFDColoring,PetscViewer);
1285: PETSC_EXTERN PetscErrorCode MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
1286: PETSC_EXTERN PetscErrorCode MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**);
1287: PETSC_EXTERN PetscErrorCode MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
1288: PETSC_EXTERN PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring);
1289: PETSC_EXTERN PetscErrorCode MatFDColoringApply(Mat,MatFDColoring,Vec,void *);
1290: PETSC_EXTERN PetscErrorCode MatFDColoringSetF(MatFDColoring,Vec);
1291: PETSC_EXTERN PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,const PetscInt*[]);
1292: PETSC_EXTERN PetscErrorCode MatFDColoringSetUp(Mat,ISColoring,MatFDColoring);
1293: PETSC_EXTERN PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring,PetscInt,PetscInt);
1294: PETSC_EXTERN PetscErrorCode MatFDColoringSetValues(Mat,MatFDColoring,const PetscScalar*);

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

1299:    Level: beginner

1301: .seealso:  MatTransposeColoringCreate()
1302: S*/
1303: typedef struct _p_MatTransposeColoring* MatTransposeColoring;

1305: PETSC_EXTERN PetscErrorCode MatTransposeColoringCreate(Mat,ISColoring,MatTransposeColoring *);
1306: PETSC_EXTERN PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring,Mat,Mat);
1307: PETSC_EXTERN PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring,Mat,Mat);
1308: PETSC_EXTERN PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring*);

1310: /*
1311:     These routines are for partitioning matrices: currently used only
1312:   for adjacency matrix, MatCreateMPIAdj().
1313: */

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

1318:    Level: beginner

1320: .seealso:  MatPartitioningCreate(), MatPartitioningType
1321: S*/
1322: typedef struct _p_MatPartitioning* MatPartitioning;

1324: /*J
1325:     MatPartitioningType - String with the name of a PETSc matrix partitioning

1327:    Level: beginner
1328: dm
1329: .seealso: MatPartitioningCreate(), MatPartitioning
1330: J*/
1331: typedef const char* MatPartitioningType;
1332: #define MATPARTITIONINGCURRENT  "current"
1333: #define MATPARTITIONINGAVERAGE   "average"
1334: #define MATPARTITIONINGSQUARE   "square"
1335: #define MATPARTITIONINGPARMETIS "parmetis"
1336: #define MATPARTITIONINGCHACO    "chaco"
1337: #define MATPARTITIONINGPARTY    "party"
1338: #define MATPARTITIONINGPTSCOTCH "ptscotch"
1339: #define MATPARTITIONINGHIERARCH  "hierarch"


1342: PETSC_EXTERN PetscErrorCode MatPartitioningCreate(MPI_Comm,MatPartitioning*);
1343: PETSC_EXTERN PetscErrorCode MatPartitioningSetType(MatPartitioning,MatPartitioningType);
1344: PETSC_EXTERN PetscErrorCode MatPartitioningSetNParts(MatPartitioning,PetscInt);
1345: PETSC_EXTERN PetscErrorCode MatPartitioningSetAdjacency(MatPartitioning,Mat);
1346: PETSC_EXTERN PetscErrorCode MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
1347: PETSC_EXTERN PetscErrorCode MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
1348: PETSC_EXTERN PetscErrorCode MatPartitioningApply(MatPartitioning,IS*);
1349: PETSC_EXTERN PetscErrorCode MatPartitioningImprove(MatPartitioning,IS*);
1350: PETSC_EXTERN PetscErrorCode MatPartitioningViewImbalance(MatPartitioning,IS);
1351: PETSC_EXTERN PetscErrorCode MatPartitioningApplyND(MatPartitioning,IS*);
1352: PETSC_EXTERN PetscErrorCode MatPartitioningDestroy(MatPartitioning*);

1354: PETSC_EXTERN PetscErrorCode MatPartitioningRegister(const char[],PetscErrorCode (*)(MatPartitioning));



1358: PETSC_EXTERN PetscErrorCode MatPartitioningView(MatPartitioning,PetscViewer);
1359: PETSC_STATIC_INLINE PetscErrorCode MatPartitioningViewFromOptions(MatPartitioning A,PetscObject B,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,B,name);}

1361: PETSC_EXTERN PetscErrorCode MatPartitioningSetFromOptions(MatPartitioning);
1362: PETSC_EXTERN PetscErrorCode MatPartitioningGetType(MatPartitioning,MatPartitioningType*);

1364: PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetRepartition(MatPartitioning part);
1365: PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
1366: PETSC_EXTERN PetscErrorCode MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *);

1368: typedef enum { MP_CHACO_MULTILEVEL=1,MP_CHACO_SPECTRAL=2,MP_CHACO_LINEAR=4,MP_CHACO_RANDOM=5,MP_CHACO_SCATTERED=6 } MPChacoGlobalType;
1369: PETSC_EXTERN const char *const MPChacoGlobalTypes[];
1370: typedef enum { MP_CHACO_KERNIGHAN=1,MP_CHACO_NONE=2 } MPChacoLocalType;
1371: PETSC_EXTERN const char *const MPChacoLocalTypes[];
1372: typedef enum { MP_CHACO_LANCZOS=0,MP_CHACO_RQI=1 } MPChacoEigenType;
1373: PETSC_EXTERN const char *const MPChacoEigenTypes[];

1375: PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetGlobal(MatPartitioning,MPChacoGlobalType);
1376: PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetGlobal(MatPartitioning,MPChacoGlobalType*);
1377: PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetLocal(MatPartitioning,MPChacoLocalType);
1378: PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetLocal(MatPartitioning,MPChacoLocalType*);
1379: PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
1380: PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
1381: PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenSolver(MatPartitioning,MPChacoEigenType*);
1382: PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenTol(MatPartitioning,PetscReal);
1383: PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenTol(MatPartitioning,PetscReal*);
1384: PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenNumber(MatPartitioning,PetscInt);
1385: PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenNumber(MatPartitioning,PetscInt*);

1387: #define MP_PARTY_OPT "opt"
1388: #define MP_PARTY_LIN "lin"
1389: #define MP_PARTY_SCA "sca"
1390: #define MP_PARTY_RAN "ran"
1391: #define MP_PARTY_GBF "gbf"
1392: #define MP_PARTY_GCF "gcf"
1393: #define MP_PARTY_BUB "bub"
1394: #define MP_PARTY_DEF "def"
1395: PETSC_EXTERN PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning,const char*);
1396: #define MP_PARTY_HELPFUL_SETS "hs"
1397: #define MP_PARTY_KERNIGHAN_LIN "kl"
1398: #define MP_PARTY_NONE "no"
1399: PETSC_EXTERN PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning,const char*);
1400: PETSC_EXTERN PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
1401: PETSC_EXTERN PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning,PetscBool);
1402: PETSC_EXTERN PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool);

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

1407: PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning,PetscReal);
1408: PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning,PetscReal*);
1409: PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning,MPPTScotchStrategyType);
1410: PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning,MPPTScotchStrategyType*);

1412: /*
1413:  * hierarchical partitioning
1414:  */
1415: PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetFineparts(MatPartitioning,IS*);
1416: PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetCoarseparts(MatPartitioning,IS*);
1417: PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNcoarseparts(MatPartitioning,PetscInt);
1418: PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNfineparts(MatPartitioning, PetscInt);

1420: PETSC_EXTERN PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*);
1421: PETSC_EXTERN PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*);

1423: /*
1424:     If you add entries here you must also add them to petsc/finclude/petscmat.h
1425: */
1426: typedef enum { MATOP_SET_VALUES=0,
1427:                MATOP_GET_ROW=1,
1428:                MATOP_RESTORE_ROW=2,
1429:                MATOP_MULT=3,
1430:                MATOP_MULT_ADD=4,
1431:                MATOP_MULT_TRANSPOSE=5,
1432:                MATOP_MULT_TRANSPOSE_ADD=6,
1433:                MATOP_SOLVE=7,
1434:                MATOP_SOLVE_ADD=8,
1435:                MATOP_SOLVE_TRANSPOSE=9,
1436:                MATOP_SOLVE_TRANSPOSE_ADD=10,
1437:                MATOP_LUFACTOR=11,
1438:                MATOP_CHOLESKYFACTOR=12,
1439:                MATOP_SOR=13,
1440:                MATOP_TRANSPOSE=14,
1441:                MATOP_GETINFO=15,
1442:                MATOP_EQUAL=16,
1443:                MATOP_GET_DIAGONAL=17,
1444:                MATOP_DIAGONAL_SCALE=18,
1445:                MATOP_NORM=19,
1446:                MATOP_ASSEMBLY_BEGIN=20,
1447:                MATOP_ASSEMBLY_END=21,
1448:                MATOP_SET_OPTION=22,
1449:                MATOP_ZERO_ENTRIES=23,
1450:                MATOP_ZERO_ROWS=24,
1451:                MATOP_LUFACTOR_SYMBOLIC=25,
1452:                MATOP_LUFACTOR_NUMERIC=26,
1453:                MATOP_CHOLESKY_FACTOR_SYMBOLIC=27,
1454:                MATOP_CHOLESKY_FACTOR_NUMERIC=28,
1455:                MATOP_SETUP_PREALLOCATION=29,
1456:                MATOP_ILUFACTOR_SYMBOLIC=30,
1457:                MATOP_ICCFACTOR_SYMBOLIC=31,
1458:                MATOP_GET_DIAGONAL_BLOCK=32,
1459:                MATOP_FREE_INTER_STRUCT=33,
1460:                MATOP_DUPLICATE=34,
1461:                MATOP_FORWARD_SOLVE=35,
1462:                MATOP_BACKWARD_SOLVE=36,
1463:                MATOP_ILUFACTOR=37,
1464:                MATOP_ICCFACTOR=38,
1465:                MATOP_AXPY=39,
1466:                MATOP_CREATE_SUBMATRICES=40,
1467:                MATOP_INCREASE_OVERLAP=41,
1468:                MATOP_GET_VALUES=42,
1469:                MATOP_COPY=43,
1470:                MATOP_GET_ROW_MAX=44,
1471:                MATOP_SCALE=45,
1472:                MATOP_SHIFT=46,
1473:                MATOP_DIAGONAL_SET=47,
1474:                MATOP_ZERO_ROWS_COLUMNS=48,
1475:                MATOP_SET_RANDOM=49,
1476:                MATOP_GET_ROW_IJ=50,
1477:                MATOP_RESTORE_ROW_IJ=51,
1478:                MATOP_GET_COLUMN_IJ=52,
1479:                MATOP_RESTORE_COLUMN_IJ=53,
1480:                MATOP_FDCOLORING_CREATE=54,
1481:                MATOP_COLORING_PATCH=55,
1482:                MATOP_SET_UNFACTORED=56,
1483:                MATOP_PERMUTE=57,
1484:                MATOP_SET_VALUES_BLOCKED=58,
1485:                MATOP_CREATE_SUBMATRIX=59,
1486:                MATOP_DESTROY=60,
1487:                MATOP_VIEW=61,
1488:                MATOP_CONVERT_FROM=62,
1489:                MATOP_MATMAT_MULT=63,
1490:                MATOP_MATMAT_MULT_SYMBOLIC=64,
1491:                MATOP_MATMAT_MULT_NUMERIC=65,
1492:                MATOP_SET_LOCAL_TO_GLOBAL_MAP=66,
1493:                MATOP_SET_VALUES_LOCAL=67,
1494:                MATOP_ZERO_ROWS_LOCAL=68,
1495:                MATOP_GET_ROW_MAX_ABS=69,
1496:                MATOP_GET_ROW_MIN_ABS=70,
1497:                MATOP_CONVERT=71,
1498:                MATOP_SET_COLORING=72,
1499:                /* MATOP_PLACEHOLDER_73=73, */
1500:                MATOP_SET_VALUES_ADIFOR=74,
1501:                MATOP_FD_COLORING_APPLY=75,
1502:                MATOP_SET_FROM_OPTIONS=76,
1503:                MATOP_MULT_CONSTRAINED=77,
1504:                MATOP_MULT_TRANSPOSE_CONSTRAIN=78,
1505:                MATOP_FIND_ZERO_DIAGONALS=79,
1506:                MATOP_MULT_MULTIPLE=80,
1507:                MATOP_SOLVE_MULTIPLE=81,
1508:                MATOP_GET_INERTIA=82,
1509:                MATOP_LOAD=83,
1510:                MATOP_IS_SYMMETRIC=84,
1511:                MATOP_IS_HERMITIAN=85,
1512:                MATOP_IS_STRUCTURALLY_SYMMETRIC=86,
1513:                MATOP_SET_VALUES_BLOCKEDLOCAL=87,
1514:                MATOP_CREATE_VECS=88,
1515:                MATOP_MAT_MULT=89,
1516:                MATOP_MAT_MULT_SYMBOLIC=90,
1517:                MATOP_MAT_MULT_NUMERIC=91,
1518:                MATOP_PTAP=92,
1519:                MATOP_PTAP_SYMBOLIC=93,
1520:                MATOP_PTAP_NUMERIC=94,
1521:                MATOP_MAT_TRANSPOSE_MULT=95,
1522:                MATOP_MAT_TRANSPOSE_MULT_SYMBO=96,
1523:                MATOP_MAT_TRANSPOSE_MULT_NUMER=97,
1524:                /* MATOP_PLACEHOLDER_98=98, */
1525:                /* MATOP_PLACEHOLDER_99=99, */
1526:                /* MATOP_PLACEHOLDER_100=100, */
1527:                /* MATOP_PLACEHOLDER_101=101, */
1528:                MATOP_CONJUGATE=102,
1529:                /* MATOP_PLACEHOLDER_103=103, */
1530:                MATOP_SET_VALUES_ROW=104,
1531:                MATOP_REAL_PART=105,
1532:                MATOP_IMAGINARY_PART=106,
1533:                MATOP_GET_ROW_UPPER_TRIANGULAR=107,
1534:                MATOP_RESTORE_ROW_UPPER_TRIANG=108,
1535:                MATOP_MAT_SOLVE=109,
1536:                MATOP_MAT_SOLVE_TRANSPOSE=110,
1537:                MATOP_GET_ROW_MIN=111,
1538:                MATOP_GET_COLUMN_VECTOR=112,
1539:                MATOP_MISSING_DIAGONAL=113,
1540:                MATOP_GET_SEQ_NONZERO_STRUCTUR=114,
1541:                MATOP_CREATE=115,
1542:                MATOP_GET_GHOSTS=116,
1543:                MATOP_GET_LOCAL_SUB_MATRIX=117,
1544:                MATOP_RESTORE_LOCALSUB_MATRIX=118,
1545:                MATOP_MULT_DIAGONAL_BLOCK=119,
1546:                MATOP_HERMITIAN_TRANSPOSE=120,
1547:                MATOP_MULT_HERMITIAN_TRANSPOSE=121,
1548:                MATOP_MULT_HERMITIAN_TRANS_ADD=122,
1549:                MATOP_GET_MULTI_PROC_BLOCK=123,
1550:                MATOP_FIND_NONZERO_ROWS=124,
1551:                MATOP_GET_COLUMN_NORMS=125,
1552:                MATOP_INVERT_BLOCK_DIAGONAL=126,
1553:                /* MATOP_PLACEHOLDER_127=127, */
1554:                MATOP_CREATE_SUB_MATRICES_MPI=128,
1555:                MATOP_SET_VALUES_BATCH=129,
1556:                MATOP_TRANSPOSE_MAT_MULT=130,
1557:                MATOP_TRANSPOSE_MAT_MULT_SYMBO=131,
1558:                MATOP_TRANSPOSE_MAT_MULT_NUMER=132,
1559:                MATOP_TRANSPOSE_COLORING_CREAT=133,
1560:                MATOP_TRANS_COLORING_APPLY_SPT=134,
1561:                MATOP_TRANS_COLORING_APPLY_DEN=135,
1562:                MATOP_RART=136,
1563:                MATOP_RART_SYMBOLIC=137,
1564:                MATOP_RART_NUMERIC=138,
1565:                MATOP_SET_BLOCK_SIZES=139,
1566:                MATOP_AYPX=140,
1567:                MATOP_RESIDUAL=141,
1568:                MATOP_FDCOLORING_SETUP=142,
1569:                MATOP_MPICONCATENATESEQ=144,
1570:                MATOP_DESTROYSUBMATRICES=145
1571:              } MatOperation;
1572: PETSC_EXTERN PetscErrorCode MatSetOperation(Mat,MatOperation,void(*)(void));
1573: PETSC_EXTERN PetscErrorCode MatGetOperation(Mat,MatOperation,void(**)(void));
1574: PETSC_EXTERN PetscErrorCode MatHasOperation(Mat,MatOperation,PetscBool *);
1575: PETSC_EXTERN PetscErrorCode MatHasCongruentLayouts(Mat,PetscBool *);
1576: PETSC_EXTERN PetscErrorCode MatFreeIntermediateDataStructures(Mat);

1578: PETSC_EXTERN PetscErrorCode MatShellSetOperation(Mat,MatOperation,void(*)(void));
1579: PETSC_EXTERN PetscErrorCode MatShellGetOperation(Mat,MatOperation,void(**)(void));
1580: PETSC_EXTERN PetscErrorCode MatShellSetContext(Mat,void*);
1581: PETSC_EXTERN PetscErrorCode MatShellTestMult(Mat,PetscErrorCode (*)(void*,Vec,Vec),Vec,void*,PetscBool*);
1582: PETSC_EXTERN PetscErrorCode MatShellTestMultTranspose(Mat,PetscErrorCode (*)(void*,Vec,Vec),Vec,void*,PetscBool*);
1583: PETSC_EXTERN PetscErrorCode MatShellSetManageScalingShifts(Mat);

1585: /*
1586:    Codes for matrices stored on disk. By default they are
1587:    stored in a universal format. By changing the format with
1588:    PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will
1589:    be stored in a way natural for the matrix, for example dense matrices
1590:    would be stored as dense. Matrices stored this way may only be
1591:    read into matrices of the same type.
1592: */
1593: #define MATRIX_BINARY_FORMAT_DENSE -1

1595: PETSC_EXTERN PetscErrorCode MatMPIBAIJSetHashTableFactor(Mat,PetscReal);

1597: PETSC_EXTERN PetscErrorCode MatISSetLocalMatType(Mat,MatType);
1598: PETSC_EXTERN PetscErrorCode MatISSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1599: PETSC_EXTERN PetscErrorCode MatISStoreL2L(Mat,PetscBool);
1600: PETSC_EXTERN PetscErrorCode MatISFixLocalEmpty(Mat,PetscBool);
1601: PETSC_EXTERN PetscErrorCode MatISGetLocalMat(Mat,Mat*);
1602: PETSC_EXTERN PetscErrorCode MatISRestoreLocalMat(Mat,Mat*);
1603: PETSC_EXTERN PetscErrorCode MatISSetLocalMat(Mat,Mat);
1604: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use the MatConvert() interface (since version 3.10)") PetscErrorCode MatISGetMPIXAIJ(Mat,MatReuse,Mat*);

1606: /*S
1607:      MatNullSpace - Object that removes a null space from a vector, i.e.
1608:          orthogonalizes the vector to a subsapce

1610:    Level: advanced

1612:   Users manual sections:
1613: .   Section 4.19 Solving Singular Systems

1615: .seealso:  MatNullSpaceCreate()
1616: S*/
1617: typedef struct _p_MatNullSpace* MatNullSpace;

1619: PETSC_EXTERN PetscErrorCode MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*);
1620: PETSC_EXTERN PetscErrorCode MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*);
1621: PETSC_EXTERN PetscErrorCode MatNullSpaceDestroy(MatNullSpace*);
1622: PETSC_EXTERN PetscErrorCode MatNullSpaceRemove(MatNullSpace,Vec);
1623: PETSC_EXTERN PetscErrorCode MatGetNullSpace(Mat, MatNullSpace *);
1624: PETSC_EXTERN PetscErrorCode MatGetTransposeNullSpace(Mat, MatNullSpace *);
1625: PETSC_EXTERN PetscErrorCode MatSetTransposeNullSpace(Mat,MatNullSpace);
1626: PETSC_EXTERN PetscErrorCode MatSetNullSpace(Mat,MatNullSpace);
1627: PETSC_EXTERN PetscErrorCode MatSetNearNullSpace(Mat,MatNullSpace);
1628: PETSC_EXTERN PetscErrorCode MatGetNearNullSpace(Mat,MatNullSpace*);
1629: PETSC_EXTERN PetscErrorCode MatNullSpaceTest(MatNullSpace,Mat,PetscBool  *);
1630: PETSC_EXTERN PetscErrorCode MatNullSpaceView(MatNullSpace,PetscViewer);
1631: PETSC_EXTERN PetscErrorCode MatNullSpaceGetVecs(MatNullSpace,PetscBool*,PetscInt*,const Vec**);
1632: PETSC_EXTERN PetscErrorCode MatNullSpaceCreateRigidBody(Vec,MatNullSpace*);

1634: PETSC_EXTERN PetscErrorCode MatReorderingSeqSBAIJ(Mat,IS);
1635: PETSC_EXTERN PetscErrorCode MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
1636: PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);
1637: PETSC_EXTERN PetscErrorCode MatSeqBAIJInvertBlockDiagonal(Mat);

1639: PETSC_EXTERN PetscErrorCode MatCreateMAIJ(Mat,PetscInt,Mat*);
1640: PETSC_EXTERN PetscErrorCode MatMAIJRedimension(Mat,PetscInt,Mat*);
1641: PETSC_EXTERN PetscErrorCode MatMAIJGetAIJ(Mat,Mat*);

1643: PETSC_EXTERN PetscErrorCode MatComputeOperator(Mat,MatType,Mat*);
1644: PETSC_EXTERN PetscErrorCode MatComputeOperatorTranspose(Mat,MatType,Mat*);

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

1649: PETSC_EXTERN PetscErrorCode MatCreateKAIJ(Mat,PetscInt,PetscInt,const PetscScalar[],const PetscScalar[],Mat*);
1650: PETSC_EXTERN PetscErrorCode MatKAIJGetAIJ(Mat,Mat*);
1651: PETSC_EXTERN PetscErrorCode MatKAIJGetS(Mat,PetscInt*,PetscInt*,PetscScalar**);
1652: PETSC_EXTERN PetscErrorCode MatKAIJGetSRead(Mat,PetscInt*,PetscInt*,const PetscScalar**);
1653: PETSC_EXTERN PetscErrorCode MatKAIJRestoreS(Mat,PetscScalar**);
1654: PETSC_EXTERN PetscErrorCode MatKAIJRestoreSRead(Mat,const PetscScalar**);
1655: PETSC_EXTERN PetscErrorCode MatKAIJGetT(Mat,PetscInt*,PetscInt*,PetscScalar**);
1656: PETSC_EXTERN PetscErrorCode MatKAIJGetTRead(Mat,PetscInt*,PetscInt*,const PetscScalar**);
1657: PETSC_EXTERN PetscErrorCode MatKAIJRestoreT(Mat,PetscScalar**);
1658: PETSC_EXTERN PetscErrorCode MatKAIJRestoreTRead(Mat,const PetscScalar**);
1659: PETSC_EXTERN PetscErrorCode MatKAIJSetAIJ(Mat,Mat);
1660: PETSC_EXTERN PetscErrorCode MatKAIJSetS(Mat,PetscInt,PetscInt,const PetscScalar []);
1661: PETSC_EXTERN PetscErrorCode MatKAIJSetT(Mat,PetscInt,PetscInt,const PetscScalar []);

1663: PETSC_EXTERN PetscErrorCode MatDiagonalScaleLocal(Mat,Vec);

1665: PETSC_EXTERN PetscErrorCode MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
1666: PETSC_EXTERN PetscErrorCode MatMFFDSetBase(Mat,Vec,Vec);
1667: PETSC_EXTERN PetscErrorCode MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*);
1668: PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*));
1669: PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec));
1670: PETSC_EXTERN PetscErrorCode MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt);
1671: PETSC_EXTERN PetscErrorCode MatMFFDResetHHistory(Mat);
1672: PETSC_EXTERN PetscErrorCode MatMFFDSetFunctionError(Mat,PetscReal);
1673: PETSC_EXTERN PetscErrorCode MatMFFDSetPeriod(Mat,PetscInt);
1674: PETSC_EXTERN PetscErrorCode MatMFFDGetH(Mat,PetscScalar *);
1675: PETSC_EXTERN PetscErrorCode MatMFFDSetOptionsPrefix(Mat,const char[]);
1676: PETSC_EXTERN PetscErrorCode MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*);
1677: PETSC_EXTERN PetscErrorCode MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*);

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

1683:     Notes:
1684:     MATMFFD is a specific MatType which uses the MatMFFD data structure

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

1688:     Level: developer

1690: .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister()
1691: S*/
1692: typedef struct _p_MatMFFD* MatMFFD;

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

1697:    Level: beginner

1699: .seealso: MatMFFDSetType(), MatMFFDRegister()
1700: J*/
1701: typedef const char* MatMFFDType;
1702: #define MATMFFD_DS  "ds"
1703: #define MATMFFD_WP  "wp"

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

1708: PETSC_EXTERN PetscErrorCode MatMFFDDSSetUmin(Mat,PetscReal);
1709: PETSC_EXTERN PetscErrorCode MatMFFDWPSetComputeNormU(Mat,PetscBool );

1711: PETSC_EXTERN PetscErrorCode MatFDColoringSetType(MatFDColoring,MatMFFDType);

1713: PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
1714: PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);

1716: /*
1717:    PETSc interface to MUMPS
1718: */
1719: #ifdef PETSC_HAVE_MUMPS
1720: PETSC_EXTERN PetscErrorCode MatMumpsSetIcntl(Mat,PetscInt,PetscInt);
1721: PETSC_EXTERN PetscErrorCode MatMumpsGetIcntl(Mat,PetscInt,PetscInt*);
1722: PETSC_EXTERN PetscErrorCode MatMumpsSetCntl(Mat,PetscInt,PetscReal);
1723: PETSC_EXTERN PetscErrorCode MatMumpsGetCntl(Mat,PetscInt,PetscReal*);

1725: PETSC_EXTERN PetscErrorCode MatMumpsGetInfo(Mat,PetscInt,PetscInt*);
1726: PETSC_EXTERN PetscErrorCode MatMumpsGetInfog(Mat,PetscInt,PetscInt*);
1727: PETSC_EXTERN PetscErrorCode MatMumpsGetRinfo(Mat,PetscInt,PetscReal*);
1728: PETSC_EXTERN PetscErrorCode MatMumpsGetRinfog(Mat,PetscInt,PetscReal*);
1729: PETSC_EXTERN PetscErrorCode MatMumpsGetInverse(Mat,Mat);
1730: PETSC_EXTERN PetscErrorCode MatMumpsGetInverseTranspose(Mat,Mat);
1731: #endif

1733: /*
1734:    PETSc interface to Mkl_Pardiso
1735: */
1736: #ifdef PETSC_HAVE_MKL_PARDISO
1737: PETSC_EXTERN PetscErrorCode MatMkl_PardisoSetCntl(Mat,PetscInt,PetscInt);
1738: #endif

1740: /*
1741:    PETSc interface to Mkl_CPardiso
1742: */
1743: #ifdef PETSC_HAVE_MKL_CPARDISO
1744: PETSC_EXTERN PetscErrorCode MatMkl_CPardisoSetCntl(Mat,PetscInt,PetscInt);
1745: #endif

1747: /*
1748:    PETSc interface to SUPERLU
1749: */
1750: #ifdef PETSC_HAVE_SUPERLU
1751: PETSC_EXTERN PetscErrorCode MatSuperluSetILUDropTol(Mat,PetscReal);
1752: #endif

1754: /*
1755:    PETSc interface to SUPERLU_DIST
1756: */
1757: #ifdef PETSC_HAVE_SUPERLU_DIST
1758: PETSC_EXTERN PetscErrorCode MatSuperluDistGetDiagU(Mat,PetscScalar*);
1759: #endif

1761: /*
1762:    PETSc interface to STRUMPACK
1763: */
1764: #ifdef PETSC_HAVE_STRUMPACK
1765: /*E
1766:     MatSTRUMPACKReordering - sparsity reducing ordering to be used in STRUMPACK

1768:     Level: intermediate
1769: E*/
1770: typedef enum {MAT_STRUMPACK_NATURAL=0,
1771:               MAT_STRUMPACK_METIS=1,
1772:               MAT_STRUMPACK_PARMETIS=2,
1773:               MAT_STRUMPACK_SCOTCH=3,
1774:               MAT_STRUMPACK_PTSCOTCH=4,
1775:               MAT_STRUMPACK_RCM=5} MatSTRUMPACKReordering;
1776: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetReordering(Mat,MatSTRUMPACKReordering);
1777: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetColPerm(Mat,PetscBool);
1778: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSRelTol(Mat,PetscReal);
1779: PETSC_DEPRECATED_FUNCTION("Use MatSTRUMPACKSetHSSRelTol() (since version 3.9)") PETSC_STATIC_INLINE PetscErrorCode MatSTRUMPACKSetHSSRelCompTol(Mat mat,PetscReal rtol) {return MatSTRUMPACKSetHSSRelTol(mat,rtol);}
1780: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSAbsTol(Mat,PetscReal);
1781: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSMinSepSize(Mat,PetscInt);
1782: PETSC_DEPRECATED_FUNCTION("Use MatSTRUMPACKSetHSSMinSepSize() (since version 3.9)") PETSC_STATIC_INLINE PetscErrorCode MatSTRUMPACKSetHSSMinSize(Mat mat,PetscInt hssminsize) {return MatSTRUMPACKSetHSSMinSepSize(mat,hssminsize);}
1783: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSMaxRank(Mat,PetscInt);
1784: PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSLeafSize(Mat,PetscInt);
1785: #endif


1788: PETSC_EXTERN PetscErrorCode MatPinToCPU(Mat,PetscBool);

1790: #ifdef PETSC_HAVE_CUDA
1791: /*E
1792:     MatCUSPARSEStorageFormat - indicates the storage format for CUSPARSE (GPU)
1793:     matrices.

1795:     Not Collective

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

1801:     Level: intermediate

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

1805: .seealso: MatCUSPARSESetFormat(), MatCUSPARSEFormatOperation
1806: E*/

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

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

1813: /*E
1814:     MatCUSPARSEFormatOperation - indicates the operation of CUSPARSE (GPU)
1815:     matrices whose operation should use a particular storage format.

1817:     Not Collective

1819: +   MAT_CUSPARSE_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult
1820: .   MAT_CUSPARSE_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult
1821: .   MAT_CUSPARSE_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult
1822: -   MAT_CUSPARSE_ALL - sets the storage format for all CUSPARSE (GPU) matrices

1824:     Level: intermediate

1826: .seealso: MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat
1827: E*/
1828: typedef enum {MAT_CUSPARSE_MULT_DIAG, MAT_CUSPARSE_MULT_OFFDIAG, MAT_CUSPARSE_MULT, MAT_CUSPARSE_ALL} MatCUSPARSEFormatOperation;

1830: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
1831: PETSC_EXTERN PetscErrorCode MatCreateAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1832: PETSC_EXTERN PetscErrorCode MatCUSPARSESetFormat(Mat,MatCUSPARSEFormatOperation,MatCUSPARSEStorageFormat);
1833: #endif

1835: #if defined(PETSC_HAVE_VIENNACL)
1836: PETSC_EXTERN PetscErrorCode MatCreateSeqAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
1837: PETSC_EXTERN PetscErrorCode MatCreateAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1838: #endif

1840: /*
1841:    PETSc interface to FFTW
1842: */
1843: #if defined(PETSC_HAVE_FFTW)
1844: PETSC_EXTERN PetscErrorCode VecScatterPetscToFFTW(Mat,Vec,Vec);
1845: PETSC_EXTERN PetscErrorCode VecScatterFFTWToPetsc(Mat,Vec,Vec);
1846: PETSC_EXTERN PetscErrorCode MatCreateVecsFFTW(Mat,Vec*,Vec*,Vec*);
1847: #endif

1849: PETSC_EXTERN PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*);
1850: PETSC_EXTERN PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*);
1851: PETSC_EXTERN PetscErrorCode MatNestGetISs(Mat,IS[],IS[]);
1852: PETSC_EXTERN PetscErrorCode MatNestGetLocalISs(Mat,IS[],IS[]);
1853: PETSC_EXTERN PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***);
1854: PETSC_EXTERN PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*);
1855: PETSC_EXTERN PetscErrorCode MatNestSetVecType(Mat,VecType);
1856: PETSC_EXTERN PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]);
1857: PETSC_EXTERN PetscErrorCode MatNestSetSubMat(Mat,PetscInt,PetscInt,Mat);

1859: PETSC_EXTERN PetscErrorCode MatChop(Mat,PetscReal);
1860: PETSC_EXTERN PetscErrorCode MatComputeBandwidth(Mat,PetscReal,PetscInt*);

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

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

1866: PETSC_INTERN PetscErrorCode MatHeaderMerge(Mat,Mat*);
1867: PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat,Mat*);

1869: #endif