Actual source code: matrix.c
petsc-3.9.4 2018-09-11
2: /*
3: This is where the abstract matrix operations are defined
4: */
6: #include <petsc/private/matimpl.h>
7: #include <petsc/private/isimpl.h>
8: #include <petsc/private/vecimpl.h>
10: /* Logging support */
11: PetscClassId MAT_CLASSID;
12: PetscClassId MAT_COLORING_CLASSID;
13: PetscClassId MAT_FDCOLORING_CLASSID;
14: PetscClassId MAT_TRANSPOSECOLORING_CLASSID;
16: PetscLogEvent MAT_Mult, MAT_Mults, MAT_MultConstrained, MAT_MultAdd, MAT_MultTranspose;
17: PetscLogEvent MAT_MultTransposeConstrained, MAT_MultTransposeAdd, MAT_Solve, MAT_Solves, MAT_SolveAdd, MAT_SolveTranspose, MAT_MatSolve;
18: PetscLogEvent MAT_SolveTransposeAdd, MAT_SOR, MAT_ForwardSolve, MAT_BackwardSolve, MAT_LUFactor, MAT_LUFactorSymbolic;
19: PetscLogEvent MAT_LUFactorNumeric, MAT_CholeskyFactor, MAT_CholeskyFactorSymbolic, MAT_CholeskyFactorNumeric, MAT_ILUFactor;
20: PetscLogEvent MAT_ILUFactorSymbolic, MAT_ICCFactorSymbolic, MAT_Copy, MAT_Convert, MAT_Scale, MAT_AssemblyBegin;
21: PetscLogEvent MAT_AssemblyEnd, MAT_SetValues, MAT_GetValues, MAT_GetRow, MAT_GetRowIJ, MAT_CreateSubMats, MAT_GetOrdering, MAT_RedundantMat, MAT_GetSeqNonzeroStructure;
22: PetscLogEvent MAT_IncreaseOverlap, MAT_Partitioning, MAT_Coarsen, MAT_ZeroEntries, MAT_Load, MAT_View, MAT_AXPY, MAT_FDColoringCreate;
23: PetscLogEvent MAT_FDColoringSetUp, MAT_FDColoringApply,MAT_Transpose,MAT_FDColoringFunction, MAT_CreateSubMat;
24: PetscLogEvent MAT_TransposeColoringCreate;
25: PetscLogEvent MAT_MatMult, MAT_MatMultSymbolic, MAT_MatMultNumeric;
26: PetscLogEvent MAT_PtAP, MAT_PtAPSymbolic, MAT_PtAPNumeric,MAT_RARt, MAT_RARtSymbolic, MAT_RARtNumeric;
27: PetscLogEvent MAT_MatTransposeMult, MAT_MatTransposeMultSymbolic, MAT_MatTransposeMultNumeric;
28: PetscLogEvent MAT_TransposeMatMult, MAT_TransposeMatMultSymbolic, MAT_TransposeMatMultNumeric;
29: PetscLogEvent MAT_MatMatMult, MAT_MatMatMultSymbolic, MAT_MatMatMultNumeric;
30: PetscLogEvent MAT_MultHermitianTranspose,MAT_MultHermitianTransposeAdd;
31: PetscLogEvent MAT_Getsymtranspose, MAT_Getsymtransreduced, MAT_Transpose_SeqAIJ, MAT_GetBrowsOfAcols;
32: PetscLogEvent MAT_GetBrowsOfAocols, MAT_Getlocalmat, MAT_Getlocalmatcondensed, MAT_Seqstompi, MAT_Seqstompinum, MAT_Seqstompisym;
33: PetscLogEvent MAT_Applypapt, MAT_Applypapt_numeric, MAT_Applypapt_symbolic, MAT_GetSequentialNonzeroStructure;
34: PetscLogEvent MAT_GetMultiProcBlock;
35: PetscLogEvent MAT_CUSPARSECopyToGPU, MAT_SetValuesBatch, MAT_SetValuesBatchI, MAT_SetValuesBatchII, MAT_SetValuesBatchIII, MAT_SetValuesBatchIV;
36: PetscLogEvent MAT_ViennaCLCopyToGPU;
37: PetscLogEvent MAT_Merge,MAT_Residual,MAT_SetRandom;
38: PetscLogEvent MATCOLORING_Apply,MATCOLORING_Comm,MATCOLORING_Local,MATCOLORING_ISCreate,MATCOLORING_SetUp,MATCOLORING_Weights;
40: const char *const MatFactorTypes[] = {"NONE","LU","CHOLESKY","ILU","ICC","ILUDT","MatFactorType","MAT_FACTOR_",0};
42: /*@
43: MatSetRandom - Sets all components of a matrix to random numbers. For sparse matrices that have been preallocated it randomly selects appropriate locations
45: Logically Collective on Mat
47: Input Parameters:
48: + x - the matrix
49: - rctx - the random number context, formed by PetscRandomCreate(), or NULL and
50: it will create one internally.
52: Output Parameter:
53: . x - the matrix
55: Example of Usage:
56: .vb
57: PetscRandomCreate(PETSC_COMM_WORLD,&rctx);
58: MatSetRandom(x,rctx);
59: PetscRandomDestroy(rctx);
60: .ve
62: Level: intermediate
64: Concepts: matrix^setting to random
65: Concepts: random^matrix
67: .seealso: MatZeroEntries(), MatSetValues(), PetscRandomCreate(), PetscRandomDestroy()
68: @*/
69: PetscErrorCode MatSetRandom(Mat x,PetscRandom rctx)
70: {
72: PetscRandom randObj = NULL;
79: if (!rctx) {
80: MPI_Comm comm;
81: PetscObjectGetComm((PetscObject)x,&comm);
82: PetscRandomCreate(comm,&randObj);
83: PetscRandomSetFromOptions(randObj);
84: rctx = randObj;
85: }
87: PetscLogEventBegin(MAT_SetRandom,x,rctx,0,0);
88: (*x->ops->setrandom)(x,rctx);
89: PetscLogEventEnd(MAT_SetRandom,x,rctx,0,0);
91: x->assembled = PETSC_TRUE;
92: PetscRandomDestroy(&randObj);
93: return(0);
94: }
96: /*@
97: MatFactorGetErrorZeroPivot - returns the pivot value that was determined to be zero and the row it occurred in
99: Logically Collective on Mat
101: Input Parameters:
102: . mat - the factored matrix
104: Output Parameter:
105: + pivot - the pivot value computed
106: - row - the row that the zero pivot occurred. Note that this row must be interpreted carefully due to row reorderings and which processes
107: the share the matrix
109: Level: advanced
111: Notes: This routine does not work for factorizations done with external packages.
112: This routine should only be called if MatGetFactorError() returns a value of MAT_FACTOR_NUMERIC_ZEROPIVOT
114: This can be called on non-factored matrices that come from, for example, matrices used in SOR.
116: .seealso: MatZeroEntries(), MatFactor(), MatGetFactor(), MatFactorSymbolic(), MatFactorClearError(), MatFactorGetErrorZeroPivot()
117: @*/
118: PetscErrorCode MatFactorGetErrorZeroPivot(Mat mat,PetscReal *pivot,PetscInt *row)
119: {
122: *pivot = mat->factorerror_zeropivot_value;
123: *row = mat->factorerror_zeropivot_row;
124: return(0);
125: }
127: /*@
128: MatFactorGetError - gets the error code from a factorization
130: Logically Collective on Mat
132: Input Parameters:
133: . mat - the factored matrix
135: Output Parameter:
136: . err - the error code
138: Level: advanced
140: Notes: This can be called on non-factored matrices that come from, for example, matrices used in SOR.
142: .seealso: MatZeroEntries(), MatFactor(), MatGetFactor(), MatFactorSymbolic(), MatFactorClearError(), MatFactorGetErrorZeroPivot()
143: @*/
144: PetscErrorCode MatFactorGetError(Mat mat,MatFactorError *err)
145: {
148: *err = mat->factorerrortype;
149: return(0);
150: }
152: /*@
153: MatFactorClearError - clears the error code in a factorization
155: Logically Collective on Mat
157: Input Parameter:
158: . mat - the factored matrix
160: Level: developer
162: Notes: This can be called on non-factored matrices that come from, for example, matrices used in SOR.
164: .seealso: MatZeroEntries(), MatFactor(), MatGetFactor(), MatFactorSymbolic(), MatFactorGetError(), MatFactorGetErrorZeroPivot()
165: @*/
166: PetscErrorCode MatFactorClearError(Mat mat)
167: {
170: mat->factorerrortype = MAT_FACTOR_NOERROR;
171: mat->factorerror_zeropivot_value = 0.0;
172: mat->factorerror_zeropivot_row = 0;
173: return(0);
174: }
176: static PetscErrorCode MatFindNonzeroRows_Basic(Mat mat,IS *keptrows)
177: {
178: PetscErrorCode ierr;
179: Vec r,l;
180: const PetscScalar *al;
181: PetscInt i,nz,gnz,N,n;
184: MatGetSize(mat,&N,NULL);
185: MatGetLocalSize(mat,&n,NULL);
186: MatCreateVecs(mat,&r,&l);
187: VecSet(l,0.0);
188: VecSetRandom(r,NULL);
189: MatMult(mat,r,l);
190: VecGetArrayRead(l,&al);
191: for (i=0,nz=0;i<n;i++) if (al[i] != 0.0) nz++;
192: MPIU_Allreduce(&nz,&gnz,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)mat));
193: if (gnz != N) {
194: PetscInt *nzr;
195: PetscMalloc1(nz,&nzr);
196: if (nz) { for (i=0,nz=0;i<n;i++) if (al[i] != 0.0) nzr[nz++] = i; }
197: ISCreateGeneral(PetscObjectComm((PetscObject)mat),nz,nzr,PETSC_OWN_POINTER,keptrows);
198: } else *keptrows = NULL;
199: VecRestoreArrayRead(l,&al);
200: VecDestroy(&l);
201: VecDestroy(&r);
202: return(0);
203: }
205: /*@
206: MatFindNonzeroRows - Locate all rows that are not completely zero in the matrix
208: Input Parameter:
209: . A - the matrix
211: Output Parameter:
212: . keptrows - the rows that are not completely zero
214: Notes: keptrows is set to NULL if all rows are nonzero.
216: Level: intermediate
218: @*/
219: PetscErrorCode MatFindNonzeroRows(Mat mat,IS *keptrows)
220: {
227: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
228: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
229: if (!mat->ops->findnonzerorows) {
230: MatFindNonzeroRows_Basic(mat,keptrows);
231: } else {
232: (*mat->ops->findnonzerorows)(mat,keptrows);
233: }
234: return(0);
235: }
237: /*@
238: MatFindZeroRows - Locate all rows that are completely zero in the matrix
240: Input Parameter:
241: . A - the matrix
243: Output Parameter:
244: . zerorows - the rows that are completely zero
246: Notes: zerorows is set to NULL if no rows are zero.
248: Level: intermediate
250: @*/
251: PetscErrorCode MatFindZeroRows(Mat mat,IS *zerorows)
252: {
254: IS keptrows;
255: PetscInt m, n;
260: MatFindNonzeroRows(mat, &keptrows);
261: /* MatFindNonzeroRows sets keptrows to NULL if there are no zero rows.
262: In keeping with this convention, we set zerorows to NULL if there are no zero
263: rows. */
264: if (keptrows == NULL) {
265: *zerorows = NULL;
266: } else {
267: MatGetOwnershipRange(mat,&m,&n);
268: ISComplement(keptrows,m,n,zerorows);
269: ISDestroy(&keptrows);
270: }
271: return(0);
272: }
274: /*@
275: MatGetDiagonalBlock - Returns the part of the matrix associated with the on-process coupling
277: Not Collective
279: Input Parameters:
280: . A - the matrix
282: Output Parameters:
283: . a - the diagonal part (which is a SEQUENTIAL matrix)
285: Notes: see the manual page for MatCreateAIJ() for more information on the "diagonal part" of the matrix.
286: Use caution, as the reference count on the returned matrix is not incremented and it is used as
287: part of the containing MPI Mat's normal operation.
289: Level: advanced
291: @*/
292: PetscErrorCode MatGetDiagonalBlock(Mat A,Mat *a)
293: {
300: if (A->factortype) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
301: if (!A->ops->getdiagonalblock) {
302: PetscMPIInt size;
303: MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);
304: if (size == 1) {
305: *a = A;
306: return(0);
307: } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Not coded for this matrix type");
308: }
309: (*A->ops->getdiagonalblock)(A,a);
310: return(0);
311: }
313: /*@
314: MatGetTrace - Gets the trace of a matrix. The sum of the diagonal entries.
316: Collective on Mat
318: Input Parameters:
319: . mat - the matrix
321: Output Parameter:
322: . trace - the sum of the diagonal entries
324: Level: advanced
326: @*/
327: PetscErrorCode MatGetTrace(Mat mat,PetscScalar *trace)
328: {
330: Vec diag;
333: MatCreateVecs(mat,&diag,NULL);
334: MatGetDiagonal(mat,diag);
335: VecSum(diag,trace);
336: VecDestroy(&diag);
337: return(0);
338: }
340: /*@
341: MatRealPart - Zeros out the imaginary part of the matrix
343: Logically Collective on Mat
345: Input Parameters:
346: . mat - the matrix
348: Level: advanced
351: .seealso: MatImaginaryPart()
352: @*/
353: PetscErrorCode MatRealPart(Mat mat)
354: {
360: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
361: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
362: if (!mat->ops->realpart) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
363: MatCheckPreallocated(mat,1);
364: (*mat->ops->realpart)(mat);
365: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_VECCUDA)
366: if (mat->valid_GPU_matrix != PETSC_OFFLOAD_UNALLOCATED) {
367: mat->valid_GPU_matrix = PETSC_OFFLOAD_CPU;
368: }
369: #endif
370: return(0);
371: }
373: /*@C
374: MatGetGhosts - Get the global index of all ghost nodes defined by the sparse matrix
376: Collective on Mat
378: Input Parameter:
379: . mat - the matrix
381: Output Parameters:
382: + nghosts - number of ghosts (note for BAIJ matrices there is one ghost for each block)
383: - ghosts - the global indices of the ghost points
385: Notes: the nghosts and ghosts are suitable to pass into VecCreateGhost()
387: Level: advanced
389: @*/
390: PetscErrorCode MatGetGhosts(Mat mat,PetscInt *nghosts,const PetscInt *ghosts[])
391: {
397: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
398: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
399: if (!mat->ops->getghosts) {
400: if (nghosts) *nghosts = 0;
401: if (ghosts) *ghosts = 0;
402: } else {
403: (*mat->ops->getghosts)(mat,nghosts,ghosts);
404: }
405: return(0);
406: }
409: /*@
410: MatImaginaryPart - Moves the imaginary part of the matrix to the real part and zeros the imaginary part
412: Logically Collective on Mat
414: Input Parameters:
415: . mat - the matrix
417: Level: advanced
420: .seealso: MatRealPart()
421: @*/
422: PetscErrorCode MatImaginaryPart(Mat mat)
423: {
429: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
430: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
431: if (!mat->ops->imaginarypart) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
432: MatCheckPreallocated(mat,1);
433: (*mat->ops->imaginarypart)(mat);
434: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_VECCUDA)
435: if (mat->valid_GPU_matrix != PETSC_OFFLOAD_UNALLOCATED) {
436: mat->valid_GPU_matrix = PETSC_OFFLOAD_CPU;
437: }
438: #endif
439: return(0);
440: }
442: /*@
443: MatMissingDiagonal - Determine if sparse matrix is missing a diagonal entry (or block entry for BAIJ matrices)
445: Not Collective
447: Input Parameter:
448: . mat - the matrix
450: Output Parameters:
451: + missing - is any diagonal missing
452: - dd - first diagonal entry that is missing (optional) on this process
454: Level: advanced
457: .seealso: MatRealPart()
458: @*/
459: PetscErrorCode MatMissingDiagonal(Mat mat,PetscBool *missing,PetscInt *dd)
460: {
466: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
467: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
468: if (!mat->ops->missingdiagonal) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
469: (*mat->ops->missingdiagonal)(mat,missing,dd);
470: return(0);
471: }
473: /*@C
474: MatGetRow - Gets a row of a matrix. You MUST call MatRestoreRow()
475: for each row that you get to ensure that your application does
476: not bleed memory.
478: Not Collective
480: Input Parameters:
481: + mat - the matrix
482: - row - the row to get
484: Output Parameters:
485: + ncols - if not NULL, the number of nonzeros in the row
486: . cols - if not NULL, the column numbers
487: - vals - if not NULL, the values
489: Notes:
490: This routine is provided for people who need to have direct access
491: to the structure of a matrix. We hope that we provide enough
492: high-level matrix routines that few users will need it.
494: MatGetRow() always returns 0-based column indices, regardless of
495: whether the internal representation is 0-based (default) or 1-based.
497: For better efficiency, set cols and/or vals to NULL if you do
498: not wish to extract these quantities.
500: The user can only examine the values extracted with MatGetRow();
501: the values cannot be altered. To change the matrix entries, one
502: must use MatSetValues().
504: You can only have one call to MatGetRow() outstanding for a particular
505: matrix at a time, per processor. MatGetRow() can only obtain rows
506: associated with the given processor, it cannot get rows from the
507: other processors; for that we suggest using MatCreateSubMatrices(), then
508: MatGetRow() on the submatrix. The row index passed to MatGetRows()
509: is in the global number of rows.
511: Fortran Notes:
512: The calling sequence from Fortran is
513: .vb
514: MatGetRow(matrix,row,ncols,cols,values,ierr)
515: Mat matrix (input)
516: integer row (input)
517: integer ncols (output)
518: integer cols(maxcols) (output)
519: double precision (or double complex) values(maxcols) output
520: .ve
521: where maxcols >= maximum nonzeros in any row of the matrix.
524: Caution:
525: Do not try to change the contents of the output arrays (cols and vals).
526: In some cases, this may corrupt the matrix.
528: Level: advanced
530: Concepts: matrices^row access
532: .seealso: MatRestoreRow(), MatSetValues(), MatGetValues(), MatCreateSubMatrices(), MatGetDiagonal()
533: @*/
534: PetscErrorCode MatGetRow(Mat mat,PetscInt row,PetscInt *ncols,const PetscInt *cols[],const PetscScalar *vals[])
535: {
537: PetscInt incols;
542: if (!mat->assembled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
543: if (mat->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
544: if (!mat->ops->getrow) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
545: MatCheckPreallocated(mat,1);
546: PetscLogEventBegin(MAT_GetRow,mat,0,0,0);
547: (*mat->ops->getrow)(mat,row,&incols,(PetscInt**)cols,(PetscScalar**)vals);
548: if (ncols) *ncols = incols;
549: PetscLogEventEnd(MAT_GetRow,mat,0,0,0);
550: return(0);
551: }
553: /*@
554: MatConjugate - replaces the matrix values with their complex conjugates
556: Logically Collective on Mat
558: Input Parameters:
559: . mat - the matrix
561: Level: advanced
563: .seealso: VecConjugate()
564: @*/
565: PetscErrorCode MatConjugate(Mat mat)
566: {
567: #if defined(PETSC_USE_COMPLEX)
572: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
573: if (!mat->ops->conjugate) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Not provided for this matrix format, send email to petsc-maint@mcs.anl.gov");
574: (*mat->ops->conjugate)(mat);
575: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_VECCUDA)
576: if (mat->valid_GPU_matrix != PETSC_OFFLOAD_UNALLOCATED) {
577: mat->valid_GPU_matrix = PETSC_OFFLOAD_CPU;
578: }
579: #endif
580: return(0);
581: #else
582: return 0;
583: #endif
584: }
586: /*@C
587: MatRestoreRow - Frees any temporary space allocated by MatGetRow().
589: Not Collective
591: Input Parameters:
592: + mat - the matrix
593: . row - the row to get
594: . ncols, cols - the number of nonzeros and their columns
595: - vals - if nonzero the column values
597: Notes:
598: This routine should be called after you have finished examining the entries.
600: This routine zeros out ncols, cols, and vals. This is to prevent accidental
601: us of the array after it has been restored. If you pass NULL, it will
602: not zero the pointers. Use of cols or vals after MatRestoreRow is invalid.
604: Fortran Notes:
605: The calling sequence from Fortran is
606: .vb
607: MatRestoreRow(matrix,row,ncols,cols,values,ierr)
608: Mat matrix (input)
609: integer row (input)
610: integer ncols (output)
611: integer cols(maxcols) (output)
612: double precision (or double complex) values(maxcols) output
613: .ve
614: Where maxcols >= maximum nonzeros in any row of the matrix.
616: In Fortran MatRestoreRow() MUST be called after MatGetRow()
617: before another call to MatGetRow() can be made.
619: Level: advanced
621: .seealso: MatGetRow()
622: @*/
623: PetscErrorCode MatRestoreRow(Mat mat,PetscInt row,PetscInt *ncols,const PetscInt *cols[],const PetscScalar *vals[])
624: {
630: if (!mat->assembled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
631: if (!mat->ops->restorerow) return(0);
632: (*mat->ops->restorerow)(mat,row,ncols,(PetscInt **)cols,(PetscScalar **)vals);
633: if (ncols) *ncols = 0;
634: if (cols) *cols = NULL;
635: if (vals) *vals = NULL;
636: return(0);
637: }
639: /*@
640: MatGetRowUpperTriangular - Sets a flag to enable calls to MatGetRow() for matrix in MATSBAIJ format.
641: You should call MatRestoreRowUpperTriangular() after calling MatGetRow/MatRestoreRow() to disable the flag.
643: Not Collective
645: Input Parameters:
646: + mat - the matrix
648: Notes:
649: The flag is to ensure that users are aware of MatGetRow() only provides the upper trianglular part of the row for the matrices in MATSBAIJ format.
651: Level: advanced
653: Concepts: matrices^row access
655: .seealso: MatRestoreRowRowUpperTriangular()
656: @*/
657: PetscErrorCode MatGetRowUpperTriangular(Mat mat)
658: {
664: if (!mat->assembled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
665: if (mat->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
666: if (!mat->ops->getrowuppertriangular) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
667: MatCheckPreallocated(mat,1);
668: (*mat->ops->getrowuppertriangular)(mat);
669: return(0);
670: }
672: /*@
673: MatRestoreRowUpperTriangular - Disable calls to MatGetRow() for matrix in MATSBAIJ format.
675: Not Collective
677: Input Parameters:
678: + mat - the matrix
680: Notes:
681: This routine should be called after you have finished MatGetRow/MatRestoreRow().
684: Level: advanced
686: .seealso: MatGetRowUpperTriangular()
687: @*/
688: PetscErrorCode MatRestoreRowUpperTriangular(Mat mat)
689: {
694: if (!mat->assembled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
695: if (!mat->ops->restorerowuppertriangular) return(0);
696: (*mat->ops->restorerowuppertriangular)(mat);
697: return(0);
698: }
700: /*@C
701: MatSetOptionsPrefix - Sets the prefix used for searching for all
702: Mat options in the database.
704: Logically Collective on Mat
706: Input Parameter:
707: + A - the Mat context
708: - prefix - the prefix to prepend to all option names
710: Notes:
711: A hyphen (-) must NOT be given at the beginning of the prefix name.
712: The first character of all runtime options is AUTOMATICALLY the hyphen.
714: Level: advanced
716: .keywords: Mat, set, options, prefix, database
718: .seealso: MatSetFromOptions()
719: @*/
720: PetscErrorCode MatSetOptionsPrefix(Mat A,const char prefix[])
721: {
726: PetscObjectSetOptionsPrefix((PetscObject)A,prefix);
727: return(0);
728: }
730: /*@C
731: MatAppendOptionsPrefix - Appends to the prefix used for searching for all
732: Mat options in the database.
734: Logically Collective on Mat
736: Input Parameters:
737: + A - the Mat context
738: - prefix - the prefix to prepend to all option names
740: Notes:
741: A hyphen (-) must NOT be given at the beginning of the prefix name.
742: The first character of all runtime options is AUTOMATICALLY the hyphen.
744: Level: advanced
746: .keywords: Mat, append, options, prefix, database
748: .seealso: MatGetOptionsPrefix()
749: @*/
750: PetscErrorCode MatAppendOptionsPrefix(Mat A,const char prefix[])
751: {
756: PetscObjectAppendOptionsPrefix((PetscObject)A,prefix);
757: return(0);
758: }
760: /*@C
761: MatGetOptionsPrefix - Sets the prefix used for searching for all
762: Mat options in the database.
764: Not Collective
766: Input Parameter:
767: . A - the Mat context
769: Output Parameter:
770: . prefix - pointer to the prefix string used
772: Notes: On the fortran side, the user should pass in a string 'prefix' of
773: sufficient length to hold the prefix.
775: Level: advanced
777: .keywords: Mat, get, options, prefix, database
779: .seealso: MatAppendOptionsPrefix()
780: @*/
781: PetscErrorCode MatGetOptionsPrefix(Mat A,const char *prefix[])
782: {
787: PetscObjectGetOptionsPrefix((PetscObject)A,prefix);
788: return(0);
789: }
791: /*@
792: MatResetPreallocation - Reset mat to use the original nonzero pattern provided by users.
794: Collective on Mat
796: Input Parameters:
797: . A - the Mat context
799: Notes:
800: The allocated memory will be shrunk after calling MatAssembly with MAT_FINAL_ASSEMBLY. Users can reset the preallocation to access the original memory.
801: Currently support MPIAIJ and SEQAIJ.
803: Level: beginner
805: .keywords: Mat, ResetPreallocation
807: .seealso: MatSeqAIJSetPreallocation(), MatMPIAIJSetPreallocation(), MatXAIJSetPreallocation()
808: @*/
809: PetscErrorCode MatResetPreallocation(Mat A)
810: {
816: PetscUseMethod(A,"MatResetPreallocation_C",(Mat),(A));
817: return(0);
818: }
821: /*@
822: MatSetUp - Sets up the internal matrix data structures for the later use.
824: Collective on Mat
826: Input Parameters:
827: . A - the Mat context
829: Notes:
830: If the user has not set preallocation for this matrix then a default preallocation that is likely to be inefficient is used.
832: If a suitable preallocation routine is used, this function does not need to be called.
834: See the Performance chapter of the PETSc users manual for how to preallocate matrices
836: Level: beginner
838: .keywords: Mat, setup
840: .seealso: MatCreate(), MatDestroy()
841: @*/
842: PetscErrorCode MatSetUp(Mat A)
843: {
844: PetscMPIInt size;
849: if (!((PetscObject)A)->type_name) {
850: MPI_Comm_size(PetscObjectComm((PetscObject)A), &size);
851: if (size == 1) {
852: MatSetType(A, MATSEQAIJ);
853: } else {
854: MatSetType(A, MATMPIAIJ);
855: }
856: }
857: if (!A->preallocated && A->ops->setup) {
858: PetscInfo(A,"Warning not preallocating matrix storage\n");
859: (*A->ops->setup)(A);
860: }
861: PetscLayoutSetUp(A->rmap);
862: PetscLayoutSetUp(A->cmap);
863: A->preallocated = PETSC_TRUE;
864: return(0);
865: }
867: #if defined(PETSC_HAVE_SAWS)
868: #include <petscviewersaws.h>
869: #endif
870: /*@C
871: MatView - Visualizes a matrix object.
873: Collective on Mat
875: Input Parameters:
876: + mat - the matrix
877: - viewer - visualization context
879: Notes:
880: The available visualization contexts include
881: + PETSC_VIEWER_STDOUT_SELF - for sequential matrices
882: . PETSC_VIEWER_STDOUT_WORLD - for parallel matrices created on PETSC_COMM_WORLD
883: . PETSC_VIEWER_STDOUT_(comm) - for matrices created on MPI communicator comm
884: - PETSC_VIEWER_DRAW_WORLD - graphical display of nonzero structure
886: The user can open alternative visualization contexts with
887: + PetscViewerASCIIOpen() - Outputs matrix to a specified file
888: . PetscViewerBinaryOpen() - Outputs matrix in binary to a
889: specified file; corresponding input uses MatLoad()
890: . PetscViewerDrawOpen() - Outputs nonzero matrix structure to
891: an X window display
892: - PetscViewerSocketOpen() - Outputs matrix to Socket viewer.
893: Currently only the sequential dense and AIJ
894: matrix types support the Socket viewer.
896: The user can call PetscViewerPushFormat() to specify the output
897: format of ASCII printed objects (when using PETSC_VIEWER_STDOUT_SELF,
898: PETSC_VIEWER_STDOUT_WORLD and PetscViewerASCIIOpen). Available formats include
899: + PETSC_VIEWER_DEFAULT - default, prints matrix contents
900: . PETSC_VIEWER_ASCII_MATLAB - prints matrix contents in Matlab format
901: . PETSC_VIEWER_ASCII_DENSE - prints entire matrix including zeros
902: . PETSC_VIEWER_ASCII_COMMON - prints matrix contents, using a sparse
903: format common among all matrix types
904: . PETSC_VIEWER_ASCII_IMPL - prints matrix contents, using an implementation-specific
905: format (which is in many cases the same as the default)
906: . PETSC_VIEWER_ASCII_INFO - prints basic information about the matrix
907: size and structure (not the matrix entries)
908: . PETSC_VIEWER_ASCII_INFO_DETAIL - prints more detailed information about
909: the matrix structure
911: Options Database Keys:
912: + -mat_view ::ascii_info - Prints info on matrix at conclusion of MatAssemblyEnd()
913: . -mat_view ::ascii_info_detail - Prints more detailed info
914: . -mat_view - Prints matrix in ASCII format
915: . -mat_view ::ascii_matlab - Prints matrix in Matlab format
916: . -mat_view draw - PetscDraws nonzero structure of matrix, using MatView() and PetscDrawOpenX().
917: . -display <name> - Sets display name (default is host)
918: . -draw_pause <sec> - Sets number of seconds to pause after display
919: . -mat_view socket - Sends matrix to socket, can be accessed from Matlab (see Users-Manual: Chapter 12 Using MATLAB with PETSc for details)
920: . -viewer_socket_machine <machine> -
921: . -viewer_socket_port <port> -
922: . -mat_view binary - save matrix to file in binary format
923: - -viewer_binary_filename <name> -
924: Level: beginner
926: Notes: see the manual page for MatLoad() for the exact format of the binary file when the binary
927: viewer is used.
929: See share/petsc/matlab/PetscBinaryRead.m for a Matlab code that can read in the binary file when the binary
930: viewer is used.
932: One can use '-mat_view draw -draw_pause -1' to pause the graphical display of matrix nonzero structure.
933: And then use the following mouse functions:
934: left mouse: zoom in
935: middle mouse: zoom out
936: right mouse: continue with the simulation
938: Concepts: matrices^viewing
939: Concepts: matrices^plotting
940: Concepts: matrices^printing
942: .seealso: PetscViewerPushFormat(), PetscViewerASCIIOpen(), PetscViewerDrawOpen(),
943: PetscViewerSocketOpen(), PetscViewerBinaryOpen(), MatLoad()
944: @*/
945: PetscErrorCode MatView(Mat mat,PetscViewer viewer)
946: {
947: PetscErrorCode ierr;
948: PetscInt rows,cols,rbs,cbs;
949: PetscBool iascii,ibinary;
950: PetscViewerFormat format;
951: PetscMPIInt size;
952: #if defined(PETSC_HAVE_SAWS)
953: PetscBool issaws;
954: #endif
959: if (!viewer) {
960: PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)mat),&viewer);
961: }
964: MatCheckPreallocated(mat,1);
965: PetscViewerGetFormat(viewer,&format);
966: MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);
967: if (size == 1 && format == PETSC_VIEWER_LOAD_BALANCE) return(0);
968: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&ibinary);
969: if (ibinary) {
970: PetscBool mpiio;
971: PetscViewerBinaryGetUseMPIIO(viewer,&mpiio);
972: if (mpiio) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"PETSc matrix viewers do not support using MPI-IO, turn off that flag");
973: }
975: PetscLogEventBegin(MAT_View,mat,viewer,0,0);
976: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
977: if ((!iascii || (format != PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL)) && mat->factortype) {
978: SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"No viewers for factored matrix except ASCII info or info_detailed");
979: }
981: #if defined(PETSC_HAVE_SAWS)
982: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSAWS,&issaws);
983: #endif
984: if (iascii) {
985: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ORDER,"Must call MatAssemblyBegin/End() before viewing matrix");
986: PetscObjectPrintClassNamePrefixType((PetscObject)mat,viewer);
987: if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
988: MatNullSpace nullsp,transnullsp;
990: PetscViewerASCIIPushTab(viewer);
991: MatGetSize(mat,&rows,&cols);
992: MatGetBlockSizes(mat,&rbs,&cbs);
993: if (rbs != 1 || cbs != 1) {
994: if (rbs != cbs) {PetscViewerASCIIPrintf(viewer,"rows=%D, cols=%D, rbs=%D, cbs = %D\n",rows,cols,rbs,cbs);}
995: else {PetscViewerASCIIPrintf(viewer,"rows=%D, cols=%D, bs=%D\n",rows,cols,rbs);}
996: } else {
997: PetscViewerASCIIPrintf(viewer,"rows=%D, cols=%D\n",rows,cols);
998: }
999: if (mat->factortype) {
1000: MatSolverType solver;
1001: MatFactorGetSolverType(mat,&solver);
1002: PetscViewerASCIIPrintf(viewer,"package used to perform factorization: %s\n",solver);
1003: }
1004: if (mat->ops->getinfo) {
1005: MatInfo info;
1006: MatGetInfo(mat,MAT_GLOBAL_SUM,&info);
1007: PetscViewerASCIIPrintf(viewer,"total: nonzeros=%.f, allocated nonzeros=%.f\n",info.nz_used,info.nz_allocated);
1008: PetscViewerASCIIPrintf(viewer,"total number of mallocs used during MatSetValues calls =%D\n",(PetscInt)info.mallocs);
1009: }
1010: MatGetNullSpace(mat,&nullsp);
1011: MatGetTransposeNullSpace(mat,&transnullsp);
1012: if (nullsp) {PetscViewerASCIIPrintf(viewer," has attached null space\n");}
1013: if (transnullsp && transnullsp != nullsp) {PetscViewerASCIIPrintf(viewer," has attached transposed null space\n");}
1014: MatGetNearNullSpace(mat,&nullsp);
1015: if (nullsp) {PetscViewerASCIIPrintf(viewer," has attached near null space\n");}
1016: }
1017: #if defined(PETSC_HAVE_SAWS)
1018: } else if (issaws) {
1019: PetscMPIInt rank;
1021: PetscObjectName((PetscObject)mat);
1022: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
1023: if (!((PetscObject)mat)->amsmem && !rank) {
1024: PetscObjectViewSAWs((PetscObject)mat,viewer);
1025: }
1026: #endif
1027: }
1028: if ((format == PETSC_VIEWER_NATIVE || format == PETSC_VIEWER_LOAD_BALANCE) && mat->ops->viewnative) {
1029: PetscViewerASCIIPushTab(viewer);
1030: (*mat->ops->viewnative)(mat,viewer);
1031: PetscViewerASCIIPopTab(viewer);
1032: } else if (mat->ops->view) {
1033: PetscViewerASCIIPushTab(viewer);
1034: (*mat->ops->view)(mat,viewer);
1035: PetscViewerASCIIPopTab(viewer);
1036: }
1037: if (iascii) {
1038: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ORDER,"Must call MatAssemblyBegin/End() before viewing matrix");
1039: PetscViewerGetFormat(viewer,&format);
1040: if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
1041: PetscViewerASCIIPopTab(viewer);
1042: }
1043: }
1044: PetscLogEventEnd(MAT_View,mat,viewer,0,0);
1045: return(0);
1046: }
1048: #if defined(PETSC_USE_DEBUG)
1049: #include <../src/sys/totalview/tv_data_display.h>
1050: PETSC_UNUSED static int TV_display_type(const struct _p_Mat *mat)
1051: {
1052: TV_add_row("Local rows", "int", &mat->rmap->n);
1053: TV_add_row("Local columns", "int", &mat->cmap->n);
1054: TV_add_row("Global rows", "int", &mat->rmap->N);
1055: TV_add_row("Global columns", "int", &mat->cmap->N);
1056: TV_add_row("Typename", TV_ascii_string_type, ((PetscObject)mat)->type_name);
1057: return TV_format_OK;
1058: }
1059: #endif
1061: /*@C
1062: MatLoad - Loads a matrix that has been stored in binary format
1063: with MatView(). The matrix format is determined from the options database.
1064: Generates a parallel MPI matrix if the communicator has more than one
1065: processor. The default matrix type is AIJ.
1067: Collective on PetscViewer
1069: Input Parameters:
1070: + newmat - the newly loaded matrix, this needs to have been created with MatCreate()
1071: or some related function before a call to MatLoad()
1072: - viewer - binary file viewer, created with PetscViewerBinaryOpen()
1074: Options Database Keys:
1075: Used with block matrix formats (MATSEQBAIJ, ...) to specify
1076: block size
1077: . -matload_block_size <bs>
1079: Level: beginner
1081: Notes:
1082: If the Mat type has not yet been given then MATAIJ is used, call MatSetFromOptions() on the
1083: Mat before calling this routine if you wish to set it from the options database.
1085: MatLoad() automatically loads into the options database any options
1086: given in the file filename.info where filename is the name of the file
1087: that was passed to the PetscViewerBinaryOpen(). The options in the info
1088: file will be ignored if you use the -viewer_binary_skip_info option.
1090: If the type or size of newmat is not set before a call to MatLoad, PETSc
1091: sets the default matrix type AIJ and sets the local and global sizes.
1092: If type and/or size is already set, then the same are used.
1094: In parallel, each processor can load a subset of rows (or the
1095: entire matrix). This routine is especially useful when a large
1096: matrix is stored on disk and only part of it is desired on each
1097: processor. For example, a parallel solver may access only some of
1098: the rows from each processor. The algorithm used here reads
1099: relatively small blocks of data rather than reading the entire
1100: matrix and then subsetting it.
1102: Notes for advanced users:
1103: Most users should not need to know the details of the binary storage
1104: format, since MatLoad() and MatView() completely hide these details.
1105: But for anyone who's interested, the standard binary matrix storage
1106: format is
1108: $ int MAT_FILE_CLASSID
1109: $ int number of rows
1110: $ int number of columns
1111: $ int total number of nonzeros
1112: $ int *number nonzeros in each row
1113: $ int *column indices of all nonzeros (starting index is zero)
1114: $ PetscScalar *values of all nonzeros
1116: PETSc automatically does the byte swapping for
1117: machines that store the bytes reversed, e.g. DEC alpha, freebsd,
1118: linux, Windows and the paragon; thus if you write your own binary
1119: read/write routines you have to swap the bytes; see PetscBinaryRead()
1120: and PetscBinaryWrite() to see how this may be done.
1122: .keywords: matrix, load, binary, input
1124: .seealso: PetscViewerBinaryOpen(), MatView(), VecLoad()
1126: @*/
1127: PetscErrorCode MatLoad(Mat newmat,PetscViewer viewer)
1128: {
1130: PetscBool isbinary,flg;
1135: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
1136: if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()");
1138: if (!((PetscObject)newmat)->type_name) {
1139: MatSetType(newmat,MATAIJ);
1140: }
1142: if (!newmat->ops->load) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatLoad is not supported for type");
1143: PetscLogEventBegin(MAT_Load,viewer,0,0,0);
1144: (*newmat->ops->load)(newmat,viewer);
1145: PetscLogEventEnd(MAT_Load,viewer,0,0,0);
1147: flg = PETSC_FALSE;
1148: PetscOptionsGetBool(((PetscObject)newmat)->options,((PetscObject)newmat)->prefix,"-matload_symmetric",&flg,NULL);
1149: if (flg) {
1150: MatSetOption(newmat,MAT_SYMMETRIC,PETSC_TRUE);
1151: MatSetOption(newmat,MAT_SYMMETRY_ETERNAL,PETSC_TRUE);
1152: }
1153: flg = PETSC_FALSE;
1154: PetscOptionsGetBool(((PetscObject)newmat)->options,((PetscObject)newmat)->prefix,"-matload_spd",&flg,NULL);
1155: if (flg) {
1156: MatSetOption(newmat,MAT_SPD,PETSC_TRUE);
1157: }
1158: return(0);
1159: }
1161: PetscErrorCode MatDestroy_Redundant(Mat_Redundant **redundant)
1162: {
1164: Mat_Redundant *redund = *redundant;
1165: PetscInt i;
1168: if (redund){
1169: if (redund->matseq) { /* via MatCreateSubMatrices() */
1170: ISDestroy(&redund->isrow);
1171: ISDestroy(&redund->iscol);
1172: MatDestroySubMatrices(1,&redund->matseq);
1173: } else {
1174: PetscFree2(redund->send_rank,redund->recv_rank);
1175: PetscFree(redund->sbuf_j);
1176: PetscFree(redund->sbuf_a);
1177: for (i=0; i<redund->nrecvs; i++) {
1178: PetscFree(redund->rbuf_j[i]);
1179: PetscFree(redund->rbuf_a[i]);
1180: }
1181: PetscFree4(redund->sbuf_nz,redund->rbuf_nz,redund->rbuf_j,redund->rbuf_a);
1182: }
1184: if (redund->subcomm) {
1185: PetscCommDestroy(&redund->subcomm);
1186: }
1187: PetscFree(redund);
1188: }
1189: return(0);
1190: }
1192: /*@
1193: MatDestroy - Frees space taken by a matrix.
1195: Collective on Mat
1197: Input Parameter:
1198: . A - the matrix
1200: Level: beginner
1202: @*/
1203: PetscErrorCode MatDestroy(Mat *A)
1204: {
1208: if (!*A) return(0);
1210: if (--((PetscObject)(*A))->refct > 0) {*A = NULL; return(0);}
1212: /* if memory was published with SAWs then destroy it */
1213: PetscObjectSAWsViewOff((PetscObject)*A);
1214: if ((*A)->ops->destroy) {
1215: (*(*A)->ops->destroy)(*A);
1216: }
1218: PetscFree((*A)->solvertype);
1219: MatDestroy_Redundant(&(*A)->redundant);
1220: MatNullSpaceDestroy(&(*A)->nullsp);
1221: MatNullSpaceDestroy(&(*A)->transnullsp);
1222: MatNullSpaceDestroy(&(*A)->nearnullsp);
1223: MatDestroy(&(*A)->schur);
1224: PetscLayoutDestroy(&(*A)->rmap);
1225: PetscLayoutDestroy(&(*A)->cmap);
1226: PetscHeaderDestroy(A);
1227: return(0);
1228: }
1230: /*@C
1231: MatSetValues - Inserts or adds a block of values into a matrix.
1232: These values may be cached, so MatAssemblyBegin() and MatAssemblyEnd()
1233: MUST be called after all calls to MatSetValues() have been completed.
1235: Not Collective
1237: Input Parameters:
1238: + mat - the matrix
1239: . v - a logically two-dimensional array of values
1240: . m, idxm - the number of rows and their global indices
1241: . n, idxn - the number of columns and their global indices
1242: - addv - either ADD_VALUES or INSERT_VALUES, where
1243: ADD_VALUES adds values to any existing entries, and
1244: INSERT_VALUES replaces existing entries with new values
1246: Notes:
1247: If you create the matrix yourself (that is not with a call to DMCreateMatrix()) then you MUST call MatXXXXSetPreallocation() or
1248: MatSetUp() before using this routine
1250: By default the values, v, are row-oriented. See MatSetOption() for other options.
1252: Calls to MatSetValues() with the INSERT_VALUES and ADD_VALUES
1253: options cannot be mixed without intervening calls to the assembly
1254: routines.
1256: MatSetValues() uses 0-based row and column numbers in Fortran
1257: as well as in C.
1259: Negative indices may be passed in idxm and idxn, these rows and columns are
1260: simply ignored. This allows easily inserting element stiffness matrices
1261: with homogeneous Dirchlet boundary conditions that you don't want represented
1262: in the matrix.
1264: Efficiency Alert:
1265: The routine MatSetValuesBlocked() may offer much better efficiency
1266: for users of block sparse formats (MATSEQBAIJ and MATMPIBAIJ).
1268: Level: beginner
1270: Developer Notes: This is labeled with C so does not automatically generate Fortran stubs and interfaces
1271: because it requires multiple Fortran interfaces depending on which arguments are scalar or arrays.
1273: Concepts: matrices^putting entries in
1275: .seealso: MatSetOption(), MatAssemblyBegin(), MatAssemblyEnd(), MatSetValuesBlocked(), MatSetValuesLocal(),
1276: InsertMode, INSERT_VALUES, ADD_VALUES
1277: @*/
1278: PetscErrorCode MatSetValues(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],const PetscScalar v[],InsertMode addv)
1279: {
1281: #if defined(PETSC_USE_DEBUG)
1282: PetscInt i,j;
1283: #endif
1288: if (!m || !n) return(0); /* no values to insert */
1292: MatCheckPreallocated(mat,1);
1293: if (mat->insertmode == NOT_SET_VALUES) {
1294: mat->insertmode = addv;
1295: }
1296: #if defined(PETSC_USE_DEBUG)
1297: else if (mat->insertmode != addv) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Cannot mix add values and insert values");
1298: if (mat->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1299: if (!mat->ops->setvalues) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
1301: for (i=0; i<m; i++) {
1302: for (j=0; j<n; j++) {
1303: if (mat->erroriffailure && PetscIsInfOrNanScalar(v[i*n+j]))
1304: #if defined(PETSC_USE_COMPLEX)
1305: SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_FP,"Inserting %g+ig at matrix entry (%D,%D)",(double)PetscRealPart(v[i*n+j]),(double)PetscImaginaryPart(v[i*n+j]),idxm[i],idxn[j]);
1306: #else
1307: SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_FP,"Inserting %g at matrix entry (%D,%D)",(double)v[i*n+j],idxm[i],idxn[j]);
1308: #endif
1309: }
1310: }
1311: #endif
1313: if (mat->assembled) {
1314: mat->was_assembled = PETSC_TRUE;
1315: mat->assembled = PETSC_FALSE;
1316: }
1317: PetscLogEventBegin(MAT_SetValues,mat,0,0,0);
1318: (*mat->ops->setvalues)(mat,m,idxm,n,idxn,v,addv);
1319: PetscLogEventEnd(MAT_SetValues,mat,0,0,0);
1320: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_VECCUDA)
1321: if (mat->valid_GPU_matrix != PETSC_OFFLOAD_UNALLOCATED) {
1322: mat->valid_GPU_matrix = PETSC_OFFLOAD_CPU;
1323: }
1324: #endif
1325: return(0);
1326: }
1329: /*@
1330: MatSetValuesRowLocal - Inserts a row (block row for BAIJ matrices) of nonzero
1331: values into a matrix
1333: Not Collective
1335: Input Parameters:
1336: + mat - the matrix
1337: . row - the (block) row to set
1338: - v - a logically two-dimensional array of values
1340: Notes:
1341: By the values, v, are column-oriented (for the block version) and sorted
1343: All the nonzeros in the row must be provided
1345: The matrix must have previously had its column indices set
1347: The row must belong to this process
1349: Level: intermediate
1351: Concepts: matrices^putting entries in
1353: .seealso: MatSetOption(), MatAssemblyBegin(), MatAssemblyEnd(), MatSetValuesBlocked(), MatSetValuesLocal(),
1354: InsertMode, INSERT_VALUES, ADD_VALUES, MatSetValues(), MatSetValuesRow(), MatSetLocalToGlobalMapping()
1355: @*/
1356: PetscErrorCode MatSetValuesRowLocal(Mat mat,PetscInt row,const PetscScalar v[])
1357: {
1359: PetscInt globalrow;
1365: ISLocalToGlobalMappingApply(mat->rmap->mapping,1,&row,&globalrow);
1366: MatSetValuesRow(mat,globalrow,v);
1367: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_VECCUDA)
1368: if (mat->valid_GPU_matrix != PETSC_OFFLOAD_UNALLOCATED) {
1369: mat->valid_GPU_matrix = PETSC_OFFLOAD_CPU;
1370: }
1371: #endif
1372: return(0);
1373: }
1375: /*@
1376: MatSetValuesRow - Inserts a row (block row for BAIJ matrices) of nonzero
1377: values into a matrix
1379: Not Collective
1381: Input Parameters:
1382: + mat - the matrix
1383: . row - the (block) row to set
1384: - v - a logically two-dimensional (column major) array of values for block matrices with blocksize larger than one, otherwise a one dimensional array of values
1386: Notes:
1387: The values, v, are column-oriented for the block version.
1389: All the nonzeros in the row must be provided
1391: THE MATRIX MUST HAVE PREVIOUSLY HAD ITS COLUMN INDICES SET. IT IS RARE THAT THIS ROUTINE IS USED, usually MatSetValues() is used.
1393: The row must belong to this process
1395: Level: advanced
1397: Concepts: matrices^putting entries in
1399: .seealso: MatSetOption(), MatAssemblyBegin(), MatAssemblyEnd(), MatSetValuesBlocked(), MatSetValuesLocal(),
1400: InsertMode, INSERT_VALUES, ADD_VALUES, MatSetValues()
1401: @*/
1402: PetscErrorCode MatSetValuesRow(Mat mat,PetscInt row,const PetscScalar v[])
1403: {
1409: MatCheckPreallocated(mat,1);
1411: #if defined(PETSC_USE_DEBUG)
1412: if (mat->insertmode == ADD_VALUES) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Cannot mix add and insert values");
1413: if (mat->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1414: #endif
1415: mat->insertmode = INSERT_VALUES;
1417: if (mat->assembled) {
1418: mat->was_assembled = PETSC_TRUE;
1419: mat->assembled = PETSC_FALSE;
1420: }
1421: PetscLogEventBegin(MAT_SetValues,mat,0,0,0);
1422: if (!mat->ops->setvaluesrow) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
1423: (*mat->ops->setvaluesrow)(mat,row,v);
1424: PetscLogEventEnd(MAT_SetValues,mat,0,0,0);
1425: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_VECCUDA)
1426: if (mat->valid_GPU_matrix != PETSC_OFFLOAD_UNALLOCATED) {
1427: mat->valid_GPU_matrix = PETSC_OFFLOAD_CPU;
1428: }
1429: #endif
1430: return(0);
1431: }
1433: /*@
1434: MatSetValuesStencil - Inserts or adds a block of values into a matrix.
1435: Using structured grid indexing
1437: Not Collective
1439: Input Parameters:
1440: + mat - the matrix
1441: . m - number of rows being entered
1442: . idxm - grid coordinates (and component number when dof > 1) for matrix rows being entered
1443: . n - number of columns being entered
1444: . idxn - grid coordinates (and component number when dof > 1) for matrix columns being entered
1445: . v - a logically two-dimensional array of values
1446: - addv - either ADD_VALUES or INSERT_VALUES, where
1447: ADD_VALUES adds values to any existing entries, and
1448: INSERT_VALUES replaces existing entries with new values
1450: Notes:
1451: By default the values, v, are row-oriented. See MatSetOption() for other options.
1453: Calls to MatSetValuesStencil() with the INSERT_VALUES and ADD_VALUES
1454: options cannot be mixed without intervening calls to the assembly
1455: routines.
1457: The grid coordinates are across the entire grid, not just the local portion
1459: MatSetValuesStencil() uses 0-based row and column numbers in Fortran
1460: as well as in C.
1462: For setting/accessing vector values via array coordinates you can use the DMDAVecGetArray() routine
1464: In order to use this routine you must either obtain the matrix with DMCreateMatrix()
1465: or call MatSetLocalToGlobalMapping() and MatSetStencil() first.
1467: The columns and rows in the stencil passed in MUST be contained within the
1468: ghost region of the given process as set with DMDACreateXXX() or MatSetStencil(). For example,
1469: if you create a DMDA with an overlap of one grid level and on a particular process its first
1470: local nonghost x logical coordinate is 6 (so its first ghost x logical coordinate is 5) the
1471: first i index you can use in your column and row indices in MatSetStencil() is 5.
1473: In Fortran idxm and idxn should be declared as
1474: $ MatStencil idxm(4,m),idxn(4,n)
1475: and the values inserted using
1476: $ idxm(MatStencil_i,1) = i
1477: $ idxm(MatStencil_j,1) = j
1478: $ idxm(MatStencil_k,1) = k
1479: $ idxm(MatStencil_c,1) = c
1480: etc
1482: For periodic boundary conditions use negative indices for values to the left (below 0; that are to be
1483: obtained by wrapping values from right edge). For values to the right of the last entry using that index plus one
1484: etc to obtain values that obtained by wrapping the values from the left edge. This does not work for anything but the
1485: DM_BOUNDARY_PERIODIC boundary type.
1487: For indices that don't mean anything for your case (like the k index when working in 2d) or the c index when you have
1488: a single value per point) you can skip filling those indices.
1490: Inspired by the structured grid interface to the HYPRE package
1491: (http://www.llnl.gov/CASC/hypre)
1493: Efficiency Alert:
1494: The routine MatSetValuesBlockedStencil() may offer much better efficiency
1495: for users of block sparse formats (MATSEQBAIJ and MATMPIBAIJ).
1497: Level: beginner
1499: Concepts: matrices^putting entries in
1501: .seealso: MatSetOption(), MatAssemblyBegin(), MatAssemblyEnd(), MatSetValuesBlocked(), MatSetValuesLocal()
1502: MatSetValues(), MatSetValuesBlockedStencil(), MatSetStencil(), DMCreateMatrix(), DMDAVecGetArray(), MatStencil
1503: @*/
1504: PetscErrorCode MatSetValuesStencil(Mat mat,PetscInt m,const MatStencil idxm[],PetscInt n,const MatStencil idxn[],const PetscScalar v[],InsertMode addv)
1505: {
1507: PetscInt buf[8192],*bufm=0,*bufn=0,*jdxm,*jdxn;
1508: PetscInt j,i,dim = mat->stencil.dim,*dims = mat->stencil.dims+1,tmp;
1509: PetscInt *starts = mat->stencil.starts,*dxm = (PetscInt*)idxm,*dxn = (PetscInt*)idxn,sdim = dim - (1 - (PetscInt)mat->stencil.noc);
1512: if (!m || !n) return(0); /* no values to insert */
1519: if ((m+n) <= (PetscInt)(sizeof(buf)/sizeof(PetscInt))) {
1520: jdxm = buf; jdxn = buf+m;
1521: } else {
1522: PetscMalloc2(m,&bufm,n,&bufn);
1523: jdxm = bufm; jdxn = bufn;
1524: }
1525: for (i=0; i<m; i++) {
1526: for (j=0; j<3-sdim; j++) dxm++;
1527: tmp = *dxm++ - starts[0];
1528: for (j=0; j<dim-1; j++) {
1529: if ((*dxm++ - starts[j+1]) < 0 || tmp < 0) tmp = -1;
1530: else tmp = tmp*dims[j] + *(dxm-1) - starts[j+1];
1531: }
1532: if (mat->stencil.noc) dxm++;
1533: jdxm[i] = tmp;
1534: }
1535: for (i=0; i<n; i++) {
1536: for (j=0; j<3-sdim; j++) dxn++;
1537: tmp = *dxn++ - starts[0];
1538: for (j=0; j<dim-1; j++) {
1539: if ((*dxn++ - starts[j+1]) < 0 || tmp < 0) tmp = -1;
1540: else tmp = tmp*dims[j] + *(dxn-1) - starts[j+1];
1541: }
1542: if (mat->stencil.noc) dxn++;
1543: jdxn[i] = tmp;
1544: }
1545: MatSetValuesLocal(mat,m,jdxm,n,jdxn,v,addv);
1546: PetscFree2(bufm,bufn);
1547: return(0);
1548: }
1550: /*@
1551: MatSetValuesBlockedStencil - Inserts or adds a block of values into a matrix.
1552: Using structured grid indexing
1554: Not Collective
1556: Input Parameters:
1557: + mat - the matrix
1558: . m - number of rows being entered
1559: . idxm - grid coordinates for matrix rows being entered
1560: . n - number of columns being entered
1561: . idxn - grid coordinates for matrix columns being entered
1562: . v - a logically two-dimensional array of values
1563: - addv - either ADD_VALUES or INSERT_VALUES, where
1564: ADD_VALUES adds values to any existing entries, and
1565: INSERT_VALUES replaces existing entries with new values
1567: Notes:
1568: By default the values, v, are row-oriented and unsorted.
1569: See MatSetOption() for other options.
1571: Calls to MatSetValuesBlockedStencil() with the INSERT_VALUES and ADD_VALUES
1572: options cannot be mixed without intervening calls to the assembly
1573: routines.
1575: The grid coordinates are across the entire grid, not just the local portion
1577: MatSetValuesBlockedStencil() uses 0-based row and column numbers in Fortran
1578: as well as in C.
1580: For setting/accessing vector values via array coordinates you can use the DMDAVecGetArray() routine
1582: In order to use this routine you must either obtain the matrix with DMCreateMatrix()
1583: or call MatSetBlockSize(), MatSetLocalToGlobalMapping() and MatSetStencil() first.
1585: The columns and rows in the stencil passed in MUST be contained within the
1586: ghost region of the given process as set with DMDACreateXXX() or MatSetStencil(). For example,
1587: if you create a DMDA with an overlap of one grid level and on a particular process its first
1588: local nonghost x logical coordinate is 6 (so its first ghost x logical coordinate is 5) the
1589: first i index you can use in your column and row indices in MatSetStencil() is 5.
1591: In Fortran idxm and idxn should be declared as
1592: $ MatStencil idxm(4,m),idxn(4,n)
1593: and the values inserted using
1594: $ idxm(MatStencil_i,1) = i
1595: $ idxm(MatStencil_j,1) = j
1596: $ idxm(MatStencil_k,1) = k
1597: etc
1599: Negative indices may be passed in idxm and idxn, these rows and columns are
1600: simply ignored. This allows easily inserting element stiffness matrices
1601: with homogeneous Dirchlet boundary conditions that you don't want represented
1602: in the matrix.
1604: Inspired by the structured grid interface to the HYPRE package
1605: (http://www.llnl.gov/CASC/hypre)
1607: Level: beginner
1609: Concepts: matrices^putting entries in
1611: .seealso: MatSetOption(), MatAssemblyBegin(), MatAssemblyEnd(), MatSetValuesBlocked(), MatSetValuesLocal()
1612: MatSetValues(), MatSetValuesStencil(), MatSetStencil(), DMCreateMatrix(), DMDAVecGetArray(), MatStencil,
1613: MatSetBlockSize(), MatSetLocalToGlobalMapping()
1614: @*/
1615: PetscErrorCode MatSetValuesBlockedStencil(Mat mat,PetscInt m,const MatStencil idxm[],PetscInt n,const MatStencil idxn[],const PetscScalar v[],InsertMode addv)
1616: {
1618: PetscInt buf[8192],*bufm=0,*bufn=0,*jdxm,*jdxn;
1619: PetscInt j,i,dim = mat->stencil.dim,*dims = mat->stencil.dims+1,tmp;
1620: PetscInt *starts = mat->stencil.starts,*dxm = (PetscInt*)idxm,*dxn = (PetscInt*)idxn,sdim = dim - (1 - (PetscInt)mat->stencil.noc);
1623: if (!m || !n) return(0); /* no values to insert */
1630: if ((m+n) <= (PetscInt)(sizeof(buf)/sizeof(PetscInt))) {
1631: jdxm = buf; jdxn = buf+m;
1632: } else {
1633: PetscMalloc2(m,&bufm,n,&bufn);
1634: jdxm = bufm; jdxn = bufn;
1635: }
1636: for (i=0; i<m; i++) {
1637: for (j=0; j<3-sdim; j++) dxm++;
1638: tmp = *dxm++ - starts[0];
1639: for (j=0; j<sdim-1; j++) {
1640: if ((*dxm++ - starts[j+1]) < 0 || tmp < 0) tmp = -1;
1641: else tmp = tmp*dims[j] + *(dxm-1) - starts[j+1];
1642: }
1643: dxm++;
1644: jdxm[i] = tmp;
1645: }
1646: for (i=0; i<n; i++) {
1647: for (j=0; j<3-sdim; j++) dxn++;
1648: tmp = *dxn++ - starts[0];
1649: for (j=0; j<sdim-1; j++) {
1650: if ((*dxn++ - starts[j+1]) < 0 || tmp < 0) tmp = -1;
1651: else tmp = tmp*dims[j] + *(dxn-1) - starts[j+1];
1652: }
1653: dxn++;
1654: jdxn[i] = tmp;
1655: }
1656: MatSetValuesBlockedLocal(mat,m,jdxm,n,jdxn,v,addv);
1657: PetscFree2(bufm,bufn);
1658: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_VECCUDA)
1659: if (mat->valid_GPU_matrix != PETSC_OFFLOAD_UNALLOCATED) {
1660: mat->valid_GPU_matrix = PETSC_OFFLOAD_CPU;
1661: }
1662: #endif
1663: return(0);
1664: }
1666: /*@
1667: MatSetStencil - Sets the grid information for setting values into a matrix via
1668: MatSetValuesStencil()
1670: Not Collective
1672: Input Parameters:
1673: + mat - the matrix
1674: . dim - dimension of the grid 1, 2, or 3
1675: . dims - number of grid points in x, y, and z direction, including ghost points on your processor
1676: . starts - starting point of ghost nodes on your processor in x, y, and z direction
1677: - dof - number of degrees of freedom per node
1680: Inspired by the structured grid interface to the HYPRE package
1681: (www.llnl.gov/CASC/hyper)
1683: For matrices generated with DMCreateMatrix() this routine is automatically called and so not needed by the
1684: user.
1686: Level: beginner
1688: Concepts: matrices^putting entries in
1690: .seealso: MatSetOption(), MatAssemblyBegin(), MatAssemblyEnd(), MatSetValuesBlocked(), MatSetValuesLocal()
1691: MatSetValues(), MatSetValuesBlockedStencil(), MatSetValuesStencil()
1692: @*/
1693: PetscErrorCode MatSetStencil(Mat mat,PetscInt dim,const PetscInt dims[],const PetscInt starts[],PetscInt dof)
1694: {
1695: PetscInt i;
1702: mat->stencil.dim = dim + (dof > 1);
1703: for (i=0; i<dim; i++) {
1704: mat->stencil.dims[i] = dims[dim-i-1]; /* copy the values in backwards */
1705: mat->stencil.starts[i] = starts[dim-i-1];
1706: }
1707: mat->stencil.dims[dim] = dof;
1708: mat->stencil.starts[dim] = 0;
1709: mat->stencil.noc = (PetscBool)(dof == 1);
1710: return(0);
1711: }
1713: /*@C
1714: MatSetValuesBlocked - Inserts or adds a block of values into a matrix.
1716: Not Collective
1718: Input Parameters:
1719: + mat - the matrix
1720: . v - a logically two-dimensional array of values
1721: . m, idxm - the number of block rows and their global block indices
1722: . n, idxn - the number of block columns and their global block indices
1723: - addv - either ADD_VALUES or INSERT_VALUES, where
1724: ADD_VALUES adds values to any existing entries, and
1725: INSERT_VALUES replaces existing entries with new values
1727: Notes:
1728: If you create the matrix yourself (that is not with a call to DMCreateMatrix()) then you MUST call
1729: MatXXXXSetPreallocation() or MatSetUp() before using this routine.
1731: The m and n count the NUMBER of blocks in the row direction and column direction,
1732: NOT the total number of rows/columns; for example, if the block size is 2 and
1733: you are passing in values for rows 2,3,4,5 then m would be 2 (not 4).
1734: The values in idxm would be 1 2; that is the first index for each block divided by
1735: the block size.
1737: Note that you must call MatSetBlockSize() when constructing this matrix (before
1738: preallocating it).
1740: By default the values, v, are row-oriented, so the layout of
1741: v is the same as for MatSetValues(). See MatSetOption() for other options.
1743: Calls to MatSetValuesBlocked() with the INSERT_VALUES and ADD_VALUES
1744: options cannot be mixed without intervening calls to the assembly
1745: routines.
1747: MatSetValuesBlocked() uses 0-based row and column numbers in Fortran
1748: as well as in C.
1750: Negative indices may be passed in idxm and idxn, these rows and columns are
1751: simply ignored. This allows easily inserting element stiffness matrices
1752: with homogeneous Dirchlet boundary conditions that you don't want represented
1753: in the matrix.
1755: Each time an entry is set within a sparse matrix via MatSetValues(),
1756: internal searching must be done to determine where to place the
1757: data in the matrix storage space. By instead inserting blocks of
1758: entries via MatSetValuesBlocked(), the overhead of matrix assembly is
1759: reduced.
1761: Example:
1762: $ Suppose m=n=2 and block size(bs) = 2 The array is
1763: $
1764: $ 1 2 | 3 4
1765: $ 5 6 | 7 8
1766: $ - - - | - - -
1767: $ 9 10 | 11 12
1768: $ 13 14 | 15 16
1769: $
1770: $ v[] should be passed in like
1771: $ v[] = [1,2,3,4,5,6,7,8,9,10,11,12,13,14,15,16]
1772: $
1773: $ If you are not using row oriented storage of v (that is you called MatSetOption(mat,MAT_ROW_ORIENTED,PETSC_FALSE)) then
1774: $ v[] = [1,5,9,13,2,6,10,14,3,7,11,15,4,8,12,16]
1776: Level: intermediate
1778: Concepts: matrices^putting entries in blocked
1780: .seealso: MatSetBlockSize(), MatSetOption(), MatAssemblyBegin(), MatAssemblyEnd(), MatSetValues(), MatSetValuesBlockedLocal()
1781: @*/
1782: PetscErrorCode MatSetValuesBlocked(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],const PetscScalar v[],InsertMode addv)
1783: {
1789: if (!m || !n) return(0); /* no values to insert */
1793: MatCheckPreallocated(mat,1);
1794: if (mat->insertmode == NOT_SET_VALUES) {
1795: mat->insertmode = addv;
1796: }
1797: #if defined(PETSC_USE_DEBUG)
1798: else if (mat->insertmode != addv) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Cannot mix add values and insert values");
1799: if (mat->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1800: if (!mat->ops->setvaluesblocked && !mat->ops->setvalues) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
1801: #endif
1803: if (mat->assembled) {
1804: mat->was_assembled = PETSC_TRUE;
1805: mat->assembled = PETSC_FALSE;
1806: }
1807: PetscLogEventBegin(MAT_SetValues,mat,0,0,0);
1808: if (mat->ops->setvaluesblocked) {
1809: (*mat->ops->setvaluesblocked)(mat,m,idxm,n,idxn,v,addv);
1810: } else {
1811: PetscInt buf[8192],*bufr=0,*bufc=0,*iidxm,*iidxn;
1812: PetscInt i,j,bs,cbs;
1813: MatGetBlockSizes(mat,&bs,&cbs);
1814: if (m*bs+n*cbs <= (PetscInt)(sizeof(buf)/sizeof(PetscInt))) {
1815: iidxm = buf; iidxn = buf + m*bs;
1816: } else {
1817: PetscMalloc2(m*bs,&bufr,n*cbs,&bufc);
1818: iidxm = bufr; iidxn = bufc;
1819: }
1820: for (i=0; i<m; i++) {
1821: for (j=0; j<bs; j++) {
1822: iidxm[i*bs+j] = bs*idxm[i] + j;
1823: }
1824: }
1825: for (i=0; i<n; i++) {
1826: for (j=0; j<cbs; j++) {
1827: iidxn[i*cbs+j] = cbs*idxn[i] + j;
1828: }
1829: }
1830: MatSetValues(mat,m*bs,iidxm,n*cbs,iidxn,v,addv);
1831: PetscFree2(bufr,bufc);
1832: }
1833: PetscLogEventEnd(MAT_SetValues,mat,0,0,0);
1834: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_VECCUDA)
1835: if (mat->valid_GPU_matrix != PETSC_OFFLOAD_UNALLOCATED) {
1836: mat->valid_GPU_matrix = PETSC_OFFLOAD_CPU;
1837: }
1838: #endif
1839: return(0);
1840: }
1842: /*@
1843: MatGetValues - Gets a block of values from a matrix.
1845: Not Collective; currently only returns a local block
1847: Input Parameters:
1848: + mat - the matrix
1849: . v - a logically two-dimensional array for storing the values
1850: . m, idxm - the number of rows and their global indices
1851: - n, idxn - the number of columns and their global indices
1853: Notes:
1854: The user must allocate space (m*n PetscScalars) for the values, v.
1855: The values, v, are then returned in a row-oriented format,
1856: analogous to that used by default in MatSetValues().
1858: MatGetValues() uses 0-based row and column numbers in
1859: Fortran as well as in C.
1861: MatGetValues() requires that the matrix has been assembled
1862: with MatAssemblyBegin()/MatAssemblyEnd(). Thus, calls to
1863: MatSetValues() and MatGetValues() CANNOT be made in succession
1864: without intermediate matrix assembly.
1866: Negative row or column indices will be ignored and those locations in v[] will be
1867: left unchanged.
1869: Level: advanced
1871: Concepts: matrices^accessing values
1873: .seealso: MatGetRow(), MatCreateSubMatrices(), MatSetValues()
1874: @*/
1875: PetscErrorCode MatGetValues(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],PetscScalar v[])
1876: {
1882: if (!m || !n) return(0);
1886: if (!mat->assembled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
1887: if (mat->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1888: if (!mat->ops->getvalues) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
1889: MatCheckPreallocated(mat,1);
1891: PetscLogEventBegin(MAT_GetValues,mat,0,0,0);
1892: (*mat->ops->getvalues)(mat,m,idxm,n,idxn,v);
1893: PetscLogEventEnd(MAT_GetValues,mat,0,0,0);
1894: return(0);
1895: }
1897: /*@
1898: MatSetValuesBatch - Adds (ADD_VALUES) many blocks of values into a matrix at once. The blocks must all be square and
1899: the same size. Currently, this can only be called once and creates the given matrix.
1901: Not Collective
1903: Input Parameters:
1904: + mat - the matrix
1905: . nb - the number of blocks
1906: . bs - the number of rows (and columns) in each block
1907: . rows - a concatenation of the rows for each block
1908: - v - a concatenation of logically two-dimensional arrays of values
1910: Notes:
1911: In the future, we will extend this routine to handle rectangular blocks, and to allow multiple calls for a given matrix.
1913: Level: advanced
1915: Concepts: matrices^putting entries in
1917: .seealso: MatSetOption(), MatAssemblyBegin(), MatAssemblyEnd(), MatSetValuesBlocked(), MatSetValuesLocal(),
1918: InsertMode, INSERT_VALUES, ADD_VALUES, MatSetValues()
1919: @*/
1920: PetscErrorCode MatSetValuesBatch(Mat mat, PetscInt nb, PetscInt bs, PetscInt rows[], const PetscScalar v[])
1921: {
1929: #if defined(PETSC_USE_DEBUG)
1930: if (mat->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1931: #endif
1933: PetscLogEventBegin(MAT_SetValuesBatch,mat,0,0,0);
1934: if (mat->ops->setvaluesbatch) {
1935: (*mat->ops->setvaluesbatch)(mat,nb,bs,rows,v);
1936: } else {
1937: PetscInt b;
1938: for (b = 0; b < nb; ++b) {
1939: MatSetValues(mat, bs, &rows[b*bs], bs, &rows[b*bs], &v[b*bs*bs], ADD_VALUES);
1940: }
1941: }
1942: PetscLogEventEnd(MAT_SetValuesBatch,mat,0,0,0);
1943: return(0);
1944: }
1946: /*@
1947: MatSetLocalToGlobalMapping - Sets a local-to-global numbering for use by
1948: the routine MatSetValuesLocal() to allow users to insert matrix entries
1949: using a local (per-processor) numbering.
1951: Not Collective
1953: Input Parameters:
1954: + x - the matrix
1955: . rmapping - row mapping created with ISLocalToGlobalMappingCreate() or ISLocalToGlobalMappingCreateIS()
1956: - cmapping - column mapping
1958: Level: intermediate
1960: Concepts: matrices^local to global mapping
1961: Concepts: local to global mapping^for matrices
1963: .seealso: MatAssemblyBegin(), MatAssemblyEnd(), MatSetValues(), MatSetValuesLocal()
1964: @*/
1965: PetscErrorCode MatSetLocalToGlobalMapping(Mat x,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping)
1966: {
1975: if (x->ops->setlocaltoglobalmapping) {
1976: (*x->ops->setlocaltoglobalmapping)(x,rmapping,cmapping);
1977: } else {
1978: PetscLayoutSetISLocalToGlobalMapping(x->rmap,rmapping);
1979: PetscLayoutSetISLocalToGlobalMapping(x->cmap,cmapping);
1980: }
1981: return(0);
1982: }
1985: /*@
1986: MatGetLocalToGlobalMapping - Gets the local-to-global numbering set by MatSetLocalToGlobalMapping()
1988: Not Collective
1990: Input Parameters:
1991: . A - the matrix
1993: Output Parameters:
1994: + rmapping - row mapping
1995: - cmapping - column mapping
1997: Level: advanced
1999: Concepts: matrices^local to global mapping
2000: Concepts: local to global mapping^for matrices
2002: .seealso: MatSetValuesLocal()
2003: @*/
2004: PetscErrorCode MatGetLocalToGlobalMapping(Mat A,ISLocalToGlobalMapping *rmapping,ISLocalToGlobalMapping *cmapping)
2005: {
2011: if (rmapping) *rmapping = A->rmap->mapping;
2012: if (cmapping) *cmapping = A->cmap->mapping;
2013: return(0);
2014: }
2016: /*@
2017: MatGetLayouts - Gets the PetscLayout objects for rows and columns
2019: Not Collective
2021: Input Parameters:
2022: . A - the matrix
2024: Output Parameters:
2025: + rmap - row layout
2026: - cmap - column layout
2028: Level: advanced
2030: .seealso: MatCreateVecs(), MatGetLocalToGlobalMapping()
2031: @*/
2032: PetscErrorCode MatGetLayouts(Mat A,PetscLayout *rmap,PetscLayout *cmap)
2033: {
2039: if (rmap) *rmap = A->rmap;
2040: if (cmap) *cmap = A->cmap;
2041: return(0);
2042: }
2044: /*@C
2045: MatSetValuesLocal - Inserts or adds values into certain locations of a matrix,
2046: using a local ordering of the nodes.
2048: Not Collective
2050: Input Parameters:
2051: + mat - the matrix
2052: . nrow, irow - number of rows and their local indices
2053: . ncol, icol - number of columns and their local indices
2054: . y - a logically two-dimensional array of values
2055: - addv - either INSERT_VALUES or ADD_VALUES, where
2056: ADD_VALUES adds values to any existing entries, and
2057: INSERT_VALUES replaces existing entries with new values
2059: Notes:
2060: If you create the matrix yourself (that is not with a call to DMCreateMatrix()) then you MUST call MatXXXXSetPreallocation() or
2061: MatSetUp() before using this routine
2063: If you create the matrix yourself (that is not with a call to DMCreateMatrix()) then you MUST call MatSetLocalToGlobalMapping() before using this routine
2065: Calls to MatSetValuesLocal() with the INSERT_VALUES and ADD_VALUES
2066: options cannot be mixed without intervening calls to the assembly
2067: routines.
2069: These values may be cached, so MatAssemblyBegin() and MatAssemblyEnd()
2070: MUST be called after all calls to MatSetValuesLocal() have been completed.
2072: Level: intermediate
2074: Concepts: matrices^putting entries in with local numbering
2076: Developer Notes: This is labeled with C so does not automatically generate Fortran stubs and interfaces
2077: because it requires multiple Fortran interfaces depending on which arguments are scalar or arrays.
2079: .seealso: MatAssemblyBegin(), MatAssemblyEnd(), MatSetValues(), MatSetLocalToGlobalMapping(),
2080: MatSetValueLocal()
2081: @*/
2082: PetscErrorCode MatSetValuesLocal(Mat mat,PetscInt nrow,const PetscInt irow[],PetscInt ncol,const PetscInt icol[],const PetscScalar y[],InsertMode addv)
2083: {
2089: MatCheckPreallocated(mat,1);
2090: if (!nrow || !ncol) return(0); /* no values to insert */
2094: if (mat->insertmode == NOT_SET_VALUES) {
2095: mat->insertmode = addv;
2096: }
2097: #if defined(PETSC_USE_DEBUG)
2098: else if (mat->insertmode != addv) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Cannot mix add values and insert values");
2099: if (mat->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2100: if (!mat->ops->setvalueslocal && !mat->ops->setvalues) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
2101: #endif
2103: if (mat->assembled) {
2104: mat->was_assembled = PETSC_TRUE;
2105: mat->assembled = PETSC_FALSE;
2106: }
2107: PetscLogEventBegin(MAT_SetValues,mat,0,0,0);
2108: if (mat->ops->setvalueslocal) {
2109: (*mat->ops->setvalueslocal)(mat,nrow,irow,ncol,icol,y,addv);
2110: } else {
2111: PetscInt buf[8192],*bufr=0,*bufc=0,*irowm,*icolm;
2112: if ((nrow+ncol) <= (PetscInt)(sizeof(buf)/sizeof(PetscInt))) {
2113: irowm = buf; icolm = buf+nrow;
2114: } else {
2115: PetscMalloc2(nrow,&bufr,ncol,&bufc);
2116: irowm = bufr; icolm = bufc;
2117: }
2118: ISLocalToGlobalMappingApply(mat->rmap->mapping,nrow,irow,irowm);
2119: ISLocalToGlobalMappingApply(mat->cmap->mapping,ncol,icol,icolm);
2120: MatSetValues(mat,nrow,irowm,ncol,icolm,y,addv);
2121: PetscFree2(bufr,bufc);
2122: }
2123: PetscLogEventEnd(MAT_SetValues,mat,0,0,0);
2124: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_VECCUDA)
2125: if (mat->valid_GPU_matrix != PETSC_OFFLOAD_UNALLOCATED) {
2126: mat->valid_GPU_matrix = PETSC_OFFLOAD_CPU;
2127: }
2128: #endif
2129: return(0);
2130: }
2132: /*@C
2133: MatSetValuesBlockedLocal - Inserts or adds values into certain locations of a matrix,
2134: using a local ordering of the nodes a block at a time.
2136: Not Collective
2138: Input Parameters:
2139: + x - the matrix
2140: . nrow, irow - number of rows and their local indices
2141: . ncol, icol - number of columns and their local indices
2142: . y - a logically two-dimensional array of values
2143: - addv - either INSERT_VALUES or ADD_VALUES, where
2144: ADD_VALUES adds values to any existing entries, and
2145: INSERT_VALUES replaces existing entries with new values
2147: Notes:
2148: If you create the matrix yourself (that is not with a call to DMCreateMatrix()) then you MUST call MatXXXXSetPreallocation() or
2149: MatSetUp() before using this routine
2151: If you create the matrix yourself (that is not with a call to DMCreateMatrix()) then you MUST call MatSetBlockSize() and MatSetLocalToGlobalMapping()
2152: before using this routineBefore calling MatSetValuesLocal(), the user must first set the
2154: Calls to MatSetValuesBlockedLocal() with the INSERT_VALUES and ADD_VALUES
2155: options cannot be mixed without intervening calls to the assembly
2156: routines.
2158: These values may be cached, so MatAssemblyBegin() and MatAssemblyEnd()
2159: MUST be called after all calls to MatSetValuesBlockedLocal() have been completed.
2161: Level: intermediate
2163: Developer Notes: This is labeled with C so does not automatically generate Fortran stubs and interfaces
2164: because it requires multiple Fortran interfaces depending on which arguments are scalar or arrays.
2166: Concepts: matrices^putting blocked values in with local numbering
2168: .seealso: MatSetBlockSize(), MatSetLocalToGlobalMapping(), MatAssemblyBegin(), MatAssemblyEnd(),
2169: MatSetValuesLocal(), MatSetValuesBlocked()
2170: @*/
2171: PetscErrorCode MatSetValuesBlockedLocal(Mat mat,PetscInt nrow,const PetscInt irow[],PetscInt ncol,const PetscInt icol[],const PetscScalar y[],InsertMode addv)
2172: {
2178: MatCheckPreallocated(mat,1);
2179: if (!nrow || !ncol) return(0); /* no values to insert */
2183: if (mat->insertmode == NOT_SET_VALUES) {
2184: mat->insertmode = addv;
2185: }
2186: #if defined(PETSC_USE_DEBUG)
2187: else if (mat->insertmode != addv) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Cannot mix add values and insert values");
2188: if (mat->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2189: if (!mat->ops->setvaluesblockedlocal && !mat->ops->setvaluesblocked && !mat->ops->setvalueslocal && !mat->ops->setvalues) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
2190: #endif
2192: if (mat->assembled) {
2193: mat->was_assembled = PETSC_TRUE;
2194: mat->assembled = PETSC_FALSE;
2195: }
2196: PetscLogEventBegin(MAT_SetValues,mat,0,0,0);
2197: if (mat->ops->setvaluesblockedlocal) {
2198: (*mat->ops->setvaluesblockedlocal)(mat,nrow,irow,ncol,icol,y,addv);
2199: } else {
2200: PetscInt buf[8192],*bufr=0,*bufc=0,*irowm,*icolm;
2201: if ((nrow+ncol) <= (PetscInt)(sizeof(buf)/sizeof(PetscInt))) {
2202: irowm = buf; icolm = buf + nrow;
2203: } else {
2204: PetscMalloc2(nrow,&bufr,ncol,&bufc);
2205: irowm = bufr; icolm = bufc;
2206: }
2207: ISLocalToGlobalMappingApplyBlock(mat->rmap->mapping,nrow,irow,irowm);
2208: ISLocalToGlobalMappingApplyBlock(mat->cmap->mapping,ncol,icol,icolm);
2209: MatSetValuesBlocked(mat,nrow,irowm,ncol,icolm,y,addv);
2210: PetscFree2(bufr,bufc);
2211: }
2212: PetscLogEventEnd(MAT_SetValues,mat,0,0,0);
2213: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_VECCUDA)
2214: if (mat->valid_GPU_matrix != PETSC_OFFLOAD_UNALLOCATED) {
2215: mat->valid_GPU_matrix = PETSC_OFFLOAD_CPU;
2216: }
2217: #endif
2218: return(0);
2219: }
2221: /*@
2222: MatMultDiagonalBlock - Computes the matrix-vector product, y = Dx. Where D is defined by the inode or block structure of the diagonal
2224: Collective on Mat and Vec
2226: Input Parameters:
2227: + mat - the matrix
2228: - x - the vector to be multiplied
2230: Output Parameters:
2231: . y - the result
2233: Notes:
2234: The vectors x and y cannot be the same. I.e., one cannot
2235: call MatMult(A,y,y).
2237: Level: developer
2239: Concepts: matrix-vector product
2241: .seealso: MatMultTranspose(), MatMultAdd(), MatMultTransposeAdd()
2242: @*/
2243: PetscErrorCode MatMultDiagonalBlock(Mat mat,Vec x,Vec y)
2244: {
2253: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2254: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2255: if (x == y) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"x and y must be different vectors");
2256: MatCheckPreallocated(mat,1);
2258: if (!mat->ops->multdiagonalblock) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"This matrix type does not have a multiply defined");
2259: (*mat->ops->multdiagonalblock)(mat,x,y);
2260: PetscObjectStateIncrease((PetscObject)y);
2261: return(0);
2262: }
2264: /* --------------------------------------------------------*/
2265: /*@
2266: MatMult - Computes the matrix-vector product, y = Ax.
2268: Neighbor-wise Collective on Mat and Vec
2270: Input Parameters:
2271: + mat - the matrix
2272: - x - the vector to be multiplied
2274: Output Parameters:
2275: . y - the result
2277: Notes:
2278: The vectors x and y cannot be the same. I.e., one cannot
2279: call MatMult(A,y,y).
2281: Level: beginner
2283: Concepts: matrix-vector product
2285: .seealso: MatMultTranspose(), MatMultAdd(), MatMultTransposeAdd()
2286: @*/
2287: PetscErrorCode MatMult(Mat mat,Vec x,Vec y)
2288: {
2296: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2297: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2298: if (x == y) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"x and y must be different vectors");
2299: #if !defined(PETSC_HAVE_CONSTRAINTS)
2300: if (mat->cmap->N != x->map->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat mat,Vec x: global dim %D %D",mat->cmap->N,x->map->N);
2301: if (mat->rmap->N != y->map->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat mat,Vec y: global dim %D %D",mat->rmap->N,y->map->N);
2302: if (mat->rmap->n != y->map->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat mat,Vec y: local dim %D %D",mat->rmap->n,y->map->n);
2303: #endif
2304: VecLocked(y,3);
2305: if (mat->erroriffailure) {VecValidValues(x,2,PETSC_TRUE);}
2306: MatCheckPreallocated(mat,1);
2308: VecLockPush(x);
2309: if (!mat->ops->mult) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"This matrix type does not have a multiply defined");
2310: PetscLogEventBegin(MAT_Mult,mat,x,y,0);
2311: (*mat->ops->mult)(mat,x,y);
2312: PetscLogEventEnd(MAT_Mult,mat,x,y,0);
2313: if (mat->erroriffailure) {VecValidValues(y,3,PETSC_FALSE);}
2314: VecLockPop(x);
2315: return(0);
2316: }
2318: /*@
2319: MatMultTranspose - Computes matrix transpose times a vector y = A^T * x.
2321: Neighbor-wise Collective on Mat and Vec
2323: Input Parameters:
2324: + mat - the matrix
2325: - x - the vector to be multiplied
2327: Output Parameters:
2328: . y - the result
2330: Notes:
2331: The vectors x and y cannot be the same. I.e., one cannot
2332: call MatMultTranspose(A,y,y).
2334: For complex numbers this does NOT compute the Hermitian (complex conjugate) transpose multiple,
2335: use MatMultHermitianTranspose()
2337: Level: beginner
2339: Concepts: matrix vector product^transpose
2341: .seealso: MatMult(), MatMultAdd(), MatMultTransposeAdd(), MatMultHermitianTranspose(), MatTranspose()
2342: @*/
2343: PetscErrorCode MatMultTranspose(Mat mat,Vec x,Vec y)
2344: {
2353: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2354: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2355: if (x == y) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"x and y must be different vectors");
2356: #if !defined(PETSC_HAVE_CONSTRAINTS)
2357: if (mat->rmap->N != x->map->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat mat,Vec x: global dim %D %D",mat->rmap->N,x->map->N);
2358: if (mat->cmap->N != y->map->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat mat,Vec y: global dim %D %D",mat->cmap->N,y->map->N);
2359: #endif
2360: if (mat->erroriffailure) {VecValidValues(x,2,PETSC_TRUE);}
2361: MatCheckPreallocated(mat,1);
2363: if (!mat->ops->multtranspose) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"This matrix type does not have a multiply transpose defined");
2364: PetscLogEventBegin(MAT_MultTranspose,mat,x,y,0);
2365: VecLockPush(x);
2366: (*mat->ops->multtranspose)(mat,x,y);
2367: VecLockPop(x);
2368: PetscLogEventEnd(MAT_MultTranspose,mat,x,y,0);
2369: PetscObjectStateIncrease((PetscObject)y);
2370: if (mat->erroriffailure) {VecValidValues(y,3,PETSC_FALSE);}
2371: return(0);
2372: }
2374: /*@
2375: MatMultHermitianTranspose - Computes matrix Hermitian transpose times a vector.
2377: Neighbor-wise Collective on Mat and Vec
2379: Input Parameters:
2380: + mat - the matrix
2381: - x - the vector to be multilplied
2383: Output Parameters:
2384: . y - the result
2386: Notes:
2387: The vectors x and y cannot be the same. I.e., one cannot
2388: call MatMultHermitianTranspose(A,y,y).
2390: Also called the conjugate transpose, complex conjugate transpose, or adjoint.
2392: For real numbers MatMultTranspose() and MatMultHermitianTranspose() are identical.
2394: Level: beginner
2396: Concepts: matrix vector product^transpose
2398: .seealso: MatMult(), MatMultAdd(), MatMultHermitianTransposeAdd(), MatMultTranspose()
2399: @*/
2400: PetscErrorCode MatMultHermitianTranspose(Mat mat,Vec x,Vec y)
2401: {
2403: Vec w;
2411: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2412: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2413: if (x == y) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"x and y must be different vectors");
2414: #if !defined(PETSC_HAVE_CONSTRAINTS)
2415: if (mat->rmap->N != x->map->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat mat,Vec x: global dim %D %D",mat->rmap->N,x->map->N);
2416: if (mat->cmap->N != y->map->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat mat,Vec y: global dim %D %D",mat->cmap->N,y->map->N);
2417: #endif
2418: MatCheckPreallocated(mat,1);
2420: PetscLogEventBegin(MAT_MultHermitianTranspose,mat,x,y,0);
2421: if (mat->ops->multhermitiantranspose) {
2422: VecLockPush(x);
2423: (*mat->ops->multhermitiantranspose)(mat,x,y);
2424: VecLockPop(x);
2425: } else {
2426: VecDuplicate(x,&w);
2427: VecCopy(x,w);
2428: VecConjugate(w);
2429: MatMultTranspose(mat,w,y);
2430: VecDestroy(&w);
2431: VecConjugate(y);
2432: }
2433: PetscLogEventEnd(MAT_MultHermitianTranspose,mat,x,y,0);
2434: PetscObjectStateIncrease((PetscObject)y);
2435: return(0);
2436: }
2438: /*@
2439: MatMultAdd - Computes v3 = v2 + A * v1.
2441: Neighbor-wise Collective on Mat and Vec
2443: Input Parameters:
2444: + mat - the matrix
2445: - v1, v2 - the vectors
2447: Output Parameters:
2448: . v3 - the result
2450: Notes:
2451: The vectors v1 and v3 cannot be the same. I.e., one cannot
2452: call MatMultAdd(A,v1,v2,v1).
2454: Level: beginner
2456: Concepts: matrix vector product^addition
2458: .seealso: MatMultTranspose(), MatMult(), MatMultTransposeAdd()
2459: @*/
2460: PetscErrorCode MatMultAdd(Mat mat,Vec v1,Vec v2,Vec v3)
2461: {
2471: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2472: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2473: if (mat->cmap->N != v1->map->N) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_SIZ,"Mat mat,Vec v1: global dim %D %D",mat->cmap->N,v1->map->N);
2474: /* if (mat->rmap->N != v2->map->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat mat,Vec v2: global dim %D %D",mat->rmap->N,v2->map->N);
2475: if (mat->rmap->N != v3->map->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat mat,Vec v3: global dim %D %D",mat->rmap->N,v3->map->N); */
2476: if (mat->rmap->n != v3->map->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat mat,Vec v3: local dim %D %D",mat->rmap->n,v3->map->n);
2477: if (mat->rmap->n != v2->map->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat mat,Vec v2: local dim %D %D",mat->rmap->n,v2->map->n);
2478: if (v1 == v3) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_IDN,"v1 and v3 must be different vectors");
2479: MatCheckPreallocated(mat,1);
2481: if (!mat->ops->multadd) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"No MatMultAdd() for matrix type '%s'",((PetscObject)mat)->type_name);
2482: PetscLogEventBegin(MAT_MultAdd,mat,v1,v2,v3);
2483: VecLockPush(v1);
2484: (*mat->ops->multadd)(mat,v1,v2,v3);
2485: VecLockPop(v1);
2486: PetscLogEventEnd(MAT_MultAdd,mat,v1,v2,v3);
2487: PetscObjectStateIncrease((PetscObject)v3);
2488: return(0);
2489: }
2491: /*@
2492: MatMultTransposeAdd - Computes v3 = v2 + A' * v1.
2494: Neighbor-wise Collective on Mat and Vec
2496: Input Parameters:
2497: + mat - the matrix
2498: - v1, v2 - the vectors
2500: Output Parameters:
2501: . v3 - the result
2503: Notes:
2504: The vectors v1 and v3 cannot be the same. I.e., one cannot
2505: call MatMultTransposeAdd(A,v1,v2,v1).
2507: Level: beginner
2509: Concepts: matrix vector product^transpose and addition
2511: .seealso: MatMultTranspose(), MatMultAdd(), MatMult()
2512: @*/
2513: PetscErrorCode MatMultTransposeAdd(Mat mat,Vec v1,Vec v2,Vec v3)
2514: {
2524: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2525: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2526: if (!mat->ops->multtransposeadd) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
2527: if (v1 == v3) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_IDN,"v1 and v3 must be different vectors");
2528: if (mat->rmap->N != v1->map->N) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_SIZ,"Mat mat,Vec v1: global dim %D %D",mat->rmap->N,v1->map->N);
2529: if (mat->cmap->N != v2->map->N) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_SIZ,"Mat mat,Vec v2: global dim %D %D",mat->cmap->N,v2->map->N);
2530: if (mat->cmap->N != v3->map->N) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_SIZ,"Mat mat,Vec v3: global dim %D %D",mat->cmap->N,v3->map->N);
2531: MatCheckPreallocated(mat,1);
2533: PetscLogEventBegin(MAT_MultTransposeAdd,mat,v1,v2,v3);
2534: VecLockPush(v1);
2535: (*mat->ops->multtransposeadd)(mat,v1,v2,v3);
2536: VecLockPop(v1);
2537: PetscLogEventEnd(MAT_MultTransposeAdd,mat,v1,v2,v3);
2538: PetscObjectStateIncrease((PetscObject)v3);
2539: return(0);
2540: }
2542: /*@
2543: MatMultHermitianTransposeAdd - Computes v3 = v2 + A^H * v1.
2545: Neighbor-wise Collective on Mat and Vec
2547: Input Parameters:
2548: + mat - the matrix
2549: - v1, v2 - the vectors
2551: Output Parameters:
2552: . v3 - the result
2554: Notes:
2555: The vectors v1 and v3 cannot be the same. I.e., one cannot
2556: call MatMultHermitianTransposeAdd(A,v1,v2,v1).
2558: Level: beginner
2560: Concepts: matrix vector product^transpose and addition
2562: .seealso: MatMultHermitianTranspose(), MatMultTranspose(), MatMultAdd(), MatMult()
2563: @*/
2564: PetscErrorCode MatMultHermitianTransposeAdd(Mat mat,Vec v1,Vec v2,Vec v3)
2565: {
2575: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2576: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2577: if (v1 == v3) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_IDN,"v1 and v3 must be different vectors");
2578: if (mat->rmap->N != v1->map->N) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_SIZ,"Mat mat,Vec v1: global dim %D %D",mat->rmap->N,v1->map->N);
2579: if (mat->cmap->N != v2->map->N) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_SIZ,"Mat mat,Vec v2: global dim %D %D",mat->cmap->N,v2->map->N);
2580: if (mat->cmap->N != v3->map->N) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_SIZ,"Mat mat,Vec v3: global dim %D %D",mat->cmap->N,v3->map->N);
2581: MatCheckPreallocated(mat,1);
2583: PetscLogEventBegin(MAT_MultHermitianTransposeAdd,mat,v1,v2,v3);
2584: VecLockPush(v1);
2585: if (mat->ops->multhermitiantransposeadd) {
2586: (*mat->ops->multhermitiantransposeadd)(mat,v1,v2,v3);
2587: } else {
2588: Vec w,z;
2589: VecDuplicate(v1,&w);
2590: VecCopy(v1,w);
2591: VecConjugate(w);
2592: VecDuplicate(v3,&z);
2593: MatMultTranspose(mat,w,z);
2594: VecDestroy(&w);
2595: VecConjugate(z);
2596: VecWAXPY(v3,1.0,v2,z);
2597: VecDestroy(&z);
2598: }
2599: VecLockPop(v1);
2600: PetscLogEventEnd(MAT_MultHermitianTransposeAdd,mat,v1,v2,v3);
2601: PetscObjectStateIncrease((PetscObject)v3);
2602: return(0);
2603: }
2605: /*@
2606: MatMultConstrained - The inner multiplication routine for a
2607: constrained matrix P^T A P.
2609: Neighbor-wise Collective on Mat and Vec
2611: Input Parameters:
2612: + mat - the matrix
2613: - x - the vector to be multilplied
2615: Output Parameters:
2616: . y - the result
2618: Notes:
2619: The vectors x and y cannot be the same. I.e., one cannot
2620: call MatMult(A,y,y).
2622: Level: beginner
2624: .keywords: matrix, multiply, matrix-vector product, constraint
2625: .seealso: MatMult(), MatMultTranspose(), MatMultAdd(), MatMultTransposeAdd()
2626: @*/
2627: PetscErrorCode MatMultConstrained(Mat mat,Vec x,Vec y)
2628: {
2635: if (!mat->assembled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2636: if (mat->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2637: if (x == y) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"x and y must be different vectors");
2638: if (mat->cmap->N != x->map->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat mat,Vec x: global dim %D %D",mat->cmap->N,x->map->N);
2639: if (mat->rmap->N != y->map->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat mat,Vec y: global dim %D %D",mat->rmap->N,y->map->N);
2640: if (mat->rmap->n != y->map->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat mat,Vec y: local dim %D %D",mat->rmap->n,y->map->n);
2642: PetscLogEventBegin(MAT_MultConstrained,mat,x,y,0);
2643: VecLockPush(x);
2644: (*mat->ops->multconstrained)(mat,x,y);
2645: VecLockPop(x);
2646: PetscLogEventEnd(MAT_MultConstrained,mat,x,y,0);
2647: PetscObjectStateIncrease((PetscObject)y);
2648: return(0);
2649: }
2651: /*@
2652: MatMultTransposeConstrained - The inner multiplication routine for a
2653: constrained matrix P^T A^T P.
2655: Neighbor-wise Collective on Mat and Vec
2657: Input Parameters:
2658: + mat - the matrix
2659: - x - the vector to be multilplied
2661: Output Parameters:
2662: . y - the result
2664: Notes:
2665: The vectors x and y cannot be the same. I.e., one cannot
2666: call MatMult(A,y,y).
2668: Level: beginner
2670: .keywords: matrix, multiply, matrix-vector product, constraint
2671: .seealso: MatMult(), MatMultTranspose(), MatMultAdd(), MatMultTransposeAdd()
2672: @*/
2673: PetscErrorCode MatMultTransposeConstrained(Mat mat,Vec x,Vec y)
2674: {
2681: if (!mat->assembled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2682: if (mat->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2683: if (x == y) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"x and y must be different vectors");
2684: if (mat->rmap->N != x->map->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat mat,Vec x: global dim %D %D",mat->cmap->N,x->map->N);
2685: if (mat->cmap->N != y->map->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat mat,Vec y: global dim %D %D",mat->rmap->N,y->map->N);
2687: PetscLogEventBegin(MAT_MultConstrained,mat,x,y,0);
2688: (*mat->ops->multtransposeconstrained)(mat,x,y);
2689: PetscLogEventEnd(MAT_MultConstrained,mat,x,y,0);
2690: PetscObjectStateIncrease((PetscObject)y);
2691: return(0);
2692: }
2694: /*@C
2695: MatGetFactorType - gets the type of factorization it is
2697: Note Collective
2698: as the flag
2700: Input Parameters:
2701: . mat - the matrix
2703: Output Parameters:
2704: . t - the type, one of MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC,MAT_FACTOR_ILUDT
2706: Level: intermediate
2708: .seealso: MatFactorType, MatGetFactor()
2709: @*/
2710: PetscErrorCode MatGetFactorType(Mat mat,MatFactorType *t)
2711: {
2715: *t = mat->factortype;
2716: return(0);
2717: }
2719: /* ------------------------------------------------------------*/
2720: /*@C
2721: MatGetInfo - Returns information about matrix storage (number of
2722: nonzeros, memory, etc.).
2724: Collective on Mat if MAT_GLOBAL_MAX or MAT_GLOBAL_SUM is used as the flag
2726: Input Parameters:
2727: . mat - the matrix
2729: Output Parameters:
2730: + flag - flag indicating the type of parameters to be returned
2731: (MAT_LOCAL - local matrix, MAT_GLOBAL_MAX - maximum over all processors,
2732: MAT_GLOBAL_SUM - sum over all processors)
2733: - info - matrix information context
2735: Notes:
2736: The MatInfo context contains a variety of matrix data, including
2737: number of nonzeros allocated and used, number of mallocs during
2738: matrix assembly, etc. Additional information for factored matrices
2739: is provided (such as the fill ratio, number of mallocs during
2740: factorization, etc.). Much of this info is printed to PETSC_STDOUT
2741: when using the runtime options
2742: $ -info -mat_view ::ascii_info
2744: Example for C/C++ Users:
2745: See the file ${PETSC_DIR}/include/petscmat.h for a complete list of
2746: data within the MatInfo context. For example,
2747: .vb
2748: MatInfo info;
2749: Mat A;
2750: double mal, nz_a, nz_u;
2752: MatGetInfo(A,MAT_LOCAL,&info);
2753: mal = info.mallocs;
2754: nz_a = info.nz_allocated;
2755: .ve
2757: Example for Fortran Users:
2758: Fortran users should declare info as a double precision
2759: array of dimension MAT_INFO_SIZE, and then extract the parameters
2760: of interest. See the file ${PETSC_DIR}/include/petsc/finclude/petscmat.h
2761: a complete list of parameter names.
2762: .vb
2763: double precision info(MAT_INFO_SIZE)
2764: double precision mal, nz_a
2765: Mat A
2766: integer ierr
2768: call MatGetInfo(A,MAT_LOCAL,info,ierr)
2769: mal = info(MAT_INFO_MALLOCS)
2770: nz_a = info(MAT_INFO_NZ_ALLOCATED)
2771: .ve
2773: Level: intermediate
2775: Concepts: matrices^getting information on
2777: Developer Note: fortran interface is not autogenerated as the f90
2778: interface defintion cannot be generated correctly [due to MatInfo]
2780: .seealso: MatStashGetInfo()
2782: @*/
2783: PetscErrorCode MatGetInfo(Mat mat,MatInfoType flag,MatInfo *info)
2784: {
2791: if (!mat->ops->getinfo) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
2792: MatCheckPreallocated(mat,1);
2793: (*mat->ops->getinfo)(mat,flag,info);
2794: return(0);
2795: }
2797: /*
2798: This is used by external packages where it is not easy to get the info from the actual
2799: matrix factorization.
2800: */
2801: PetscErrorCode MatGetInfo_External(Mat A,MatInfoType flag,MatInfo *info)
2802: {
2806: PetscMemzero(info,sizeof(MatInfo));
2807: return(0);
2808: }
2810: /* ----------------------------------------------------------*/
2812: /*@C
2813: MatLUFactor - Performs in-place LU factorization of matrix.
2815: Collective on Mat
2817: Input Parameters:
2818: + mat - the matrix
2819: . row - row permutation
2820: . col - column permutation
2821: - info - options for factorization, includes
2822: $ fill - expected fill as ratio of original fill.
2823: $ dtcol - pivot tolerance (0 no pivot, 1 full column pivoting)
2824: $ Run with the option -info to determine an optimal value to use
2826: Notes:
2827: Most users should employ the simplified KSP interface for linear solvers
2828: instead of working directly with matrix algebra routines such as this.
2829: See, e.g., KSPCreate().
2831: This changes the state of the matrix to a factored matrix; it cannot be used
2832: for example with MatSetValues() unless one first calls MatSetUnfactored().
2834: Level: developer
2836: Concepts: matrices^LU factorization
2838: .seealso: MatLUFactorSymbolic(), MatLUFactorNumeric(), MatCholeskyFactor(),
2839: MatGetOrdering(), MatSetUnfactored(), MatFactorInfo, MatGetFactor()
2841: Developer Note: fortran interface is not autogenerated as the f90
2842: interface defintion cannot be generated correctly [due to MatFactorInfo]
2844: @*/
2845: PetscErrorCode MatLUFactor(Mat mat,IS row,IS col,const MatFactorInfo *info)
2846: {
2848: MatFactorInfo tinfo;
2856: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2857: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2858: if (!mat->ops->lufactor) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
2859: MatCheckPreallocated(mat,1);
2860: if (!info) {
2861: MatFactorInfoInitialize(&tinfo);
2862: info = &tinfo;
2863: }
2865: PetscLogEventBegin(MAT_LUFactor,mat,row,col,0);
2866: (*mat->ops->lufactor)(mat,row,col,info);
2867: PetscLogEventEnd(MAT_LUFactor,mat,row,col,0);
2868: PetscObjectStateIncrease((PetscObject)mat);
2869: return(0);
2870: }
2872: /*@C
2873: MatILUFactor - Performs in-place ILU factorization of matrix.
2875: Collective on Mat
2877: Input Parameters:
2878: + mat - the matrix
2879: . row - row permutation
2880: . col - column permutation
2881: - info - structure containing
2882: $ levels - number of levels of fill.
2883: $ expected fill - as ratio of original fill.
2884: $ 1 or 0 - indicating force fill on diagonal (improves robustness for matrices
2885: missing diagonal entries)
2887: Notes:
2888: Probably really in-place only when level of fill is zero, otherwise allocates
2889: new space to store factored matrix and deletes previous memory.
2891: Most users should employ the simplified KSP interface for linear solvers
2892: instead of working directly with matrix algebra routines such as this.
2893: See, e.g., KSPCreate().
2895: Level: developer
2897: Concepts: matrices^ILU factorization
2899: .seealso: MatILUFactorSymbolic(), MatLUFactorNumeric(), MatCholeskyFactor(), MatFactorInfo
2901: Developer Note: fortran interface is not autogenerated as the f90
2902: interface defintion cannot be generated correctly [due to MatFactorInfo]
2904: @*/
2905: PetscErrorCode MatILUFactor(Mat mat,IS row,IS col,const MatFactorInfo *info)
2906: {
2915: if (mat->rmap->N != mat->cmap->N) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONG,"matrix must be square");
2916: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2917: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2918: if (!mat->ops->ilufactor) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
2919: MatCheckPreallocated(mat,1);
2921: PetscLogEventBegin(MAT_ILUFactor,mat,row,col,0);
2922: (*mat->ops->ilufactor)(mat,row,col,info);
2923: PetscLogEventEnd(MAT_ILUFactor,mat,row,col,0);
2924: PetscObjectStateIncrease((PetscObject)mat);
2925: return(0);
2926: }
2928: /*@C
2929: MatLUFactorSymbolic - Performs symbolic LU factorization of matrix.
2930: Call this routine before calling MatLUFactorNumeric().
2932: Collective on Mat
2934: Input Parameters:
2935: + fact - the factor matrix obtained with MatGetFactor()
2936: . mat - the matrix
2937: . row, col - row and column permutations
2938: - info - options for factorization, includes
2939: $ fill - expected fill as ratio of original fill.
2940: $ dtcol - pivot tolerance (0 no pivot, 1 full column pivoting)
2941: $ Run with the option -info to determine an optimal value to use
2944: Notes: See Users-Manual: ch_mat for additional information about choosing the fill factor for better efficiency.
2946: Most users should employ the simplified KSP interface for linear solvers
2947: instead of working directly with matrix algebra routines such as this.
2948: See, e.g., KSPCreate().
2950: Level: developer
2952: Concepts: matrices^LU symbolic factorization
2954: .seealso: MatLUFactor(), MatLUFactorNumeric(), MatCholeskyFactor(), MatFactorInfo, MatFactorInfoInitialize()
2956: Developer Note: fortran interface is not autogenerated as the f90
2957: interface defintion cannot be generated correctly [due to MatFactorInfo]
2959: @*/
2960: PetscErrorCode MatLUFactorSymbolic(Mat fact,Mat mat,IS row,IS col,const MatFactorInfo *info)
2961: {
2971: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2972: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2973: if (!(fact)->ops->lufactorsymbolic) {
2974: MatSolverType spackage;
2975: MatFactorGetSolverType(fact,&spackage);
2976: SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Matrix type %s symbolic LU using solver package %s",((PetscObject)mat)->type_name,spackage);
2977: }
2978: MatCheckPreallocated(mat,2);
2980: PetscLogEventBegin(MAT_LUFactorSymbolic,mat,row,col,0);
2981: (fact->ops->lufactorsymbolic)(fact,mat,row,col,info);
2982: PetscLogEventEnd(MAT_LUFactorSymbolic,mat,row,col,0);
2983: PetscObjectStateIncrease((PetscObject)fact);
2984: return(0);
2985: }
2987: /*@C
2988: MatLUFactorNumeric - Performs numeric LU factorization of a matrix.
2989: Call this routine after first calling MatLUFactorSymbolic().
2991: Collective on Mat
2993: Input Parameters:
2994: + fact - the factor matrix obtained with MatGetFactor()
2995: . mat - the matrix
2996: - info - options for factorization
2998: Notes:
2999: See MatLUFactor() for in-place factorization. See
3000: MatCholeskyFactorNumeric() for the symmetric, positive definite case.
3002: Most users should employ the simplified KSP interface for linear solvers
3003: instead of working directly with matrix algebra routines such as this.
3004: See, e.g., KSPCreate().
3006: Level: developer
3008: Concepts: matrices^LU numeric factorization
3010: .seealso: MatLUFactorSymbolic(), MatLUFactor(), MatCholeskyFactor()
3012: Developer Note: fortran interface is not autogenerated as the f90
3013: interface defintion cannot be generated correctly [due to MatFactorInfo]
3015: @*/
3016: PetscErrorCode MatLUFactorNumeric(Mat fact,Mat mat,const MatFactorInfo *info)
3017: {
3025: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
3026: if (mat->rmap->N != (fact)->rmap->N || mat->cmap->N != (fact)->cmap->N) SETERRQ4(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_SIZ,"Mat mat,Mat fact: global dimensions are different %D should = %D %D should = %D",mat->rmap->N,(fact)->rmap->N,mat->cmap->N,(fact)->cmap->N);
3028: if (!(fact)->ops->lufactornumeric) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s numeric LU",((PetscObject)mat)->type_name);
3029: MatCheckPreallocated(mat,2);
3030: PetscLogEventBegin(MAT_LUFactorNumeric,mat,fact,0,0);
3031: (fact->ops->lufactornumeric)(fact,mat,info);
3032: PetscLogEventEnd(MAT_LUFactorNumeric,mat,fact,0,0);
3033: MatViewFromOptions(fact,NULL,"-mat_factor_view");
3034: PetscObjectStateIncrease((PetscObject)fact);
3035: return(0);
3036: }
3038: /*@C
3039: MatCholeskyFactor - Performs in-place Cholesky factorization of a
3040: symmetric matrix.
3042: Collective on Mat
3044: Input Parameters:
3045: + mat - the matrix
3046: . perm - row and column permutations
3047: - f - expected fill as ratio of original fill
3049: Notes:
3050: See MatLUFactor() for the nonsymmetric case. See also
3051: MatCholeskyFactorSymbolic(), and MatCholeskyFactorNumeric().
3053: Most users should employ the simplified KSP interface for linear solvers
3054: instead of working directly with matrix algebra routines such as this.
3055: See, e.g., KSPCreate().
3057: Level: developer
3059: Concepts: matrices^Cholesky factorization
3061: .seealso: MatLUFactor(), MatCholeskyFactorSymbolic(), MatCholeskyFactorNumeric()
3062: MatGetOrdering()
3064: Developer Note: fortran interface is not autogenerated as the f90
3065: interface defintion cannot be generated correctly [due to MatFactorInfo]
3067: @*/
3068: PetscErrorCode MatCholeskyFactor(Mat mat,IS perm,const MatFactorInfo *info)
3069: {
3077: if (mat->rmap->N != mat->cmap->N) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONG,"Matrix must be square");
3078: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
3079: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
3080: if (!mat->ops->choleskyfactor) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"In-place factorization for Mat type %s is not supported, try out-of-place factorization. See MatCholeskyFactorSymbolic/Numeric",((PetscObject)mat)->type_name);
3081: MatCheckPreallocated(mat,1);
3083: PetscLogEventBegin(MAT_CholeskyFactor,mat,perm,0,0);
3084: (*mat->ops->choleskyfactor)(mat,perm,info);
3085: PetscLogEventEnd(MAT_CholeskyFactor,mat,perm,0,0);
3086: PetscObjectStateIncrease((PetscObject)mat);
3087: return(0);
3088: }
3090: /*@C
3091: MatCholeskyFactorSymbolic - Performs symbolic Cholesky factorization
3092: of a symmetric matrix.
3094: Collective on Mat
3096: Input Parameters:
3097: + fact - the factor matrix obtained with MatGetFactor()
3098: . mat - the matrix
3099: . perm - row and column permutations
3100: - info - options for factorization, includes
3101: $ fill - expected fill as ratio of original fill.
3102: $ dtcol - pivot tolerance (0 no pivot, 1 full column pivoting)
3103: $ Run with the option -info to determine an optimal value to use
3105: Notes:
3106: See MatLUFactorSymbolic() for the nonsymmetric case. See also
3107: MatCholeskyFactor() and MatCholeskyFactorNumeric().
3109: Most users should employ the simplified KSP interface for linear solvers
3110: instead of working directly with matrix algebra routines such as this.
3111: See, e.g., KSPCreate().
3113: Level: developer
3115: Concepts: matrices^Cholesky symbolic factorization
3117: .seealso: MatLUFactorSymbolic(), MatCholeskyFactor(), MatCholeskyFactorNumeric()
3118: MatGetOrdering()
3120: Developer Note: fortran interface is not autogenerated as the f90
3121: interface defintion cannot be generated correctly [due to MatFactorInfo]
3123: @*/
3124: PetscErrorCode MatCholeskyFactorSymbolic(Mat fact,Mat mat,IS perm,const MatFactorInfo *info)
3125: {
3134: if (mat->rmap->N != mat->cmap->N) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONG,"Matrix must be square");
3135: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
3136: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
3137: if (!(fact)->ops->choleskyfactorsymbolic) {
3138: MatSolverType spackage;
3139: MatFactorGetSolverType(fact,&spackage);
3140: SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s symbolic factor Cholesky using solver package %s",((PetscObject)mat)->type_name,spackage);
3141: }
3142: MatCheckPreallocated(mat,2);
3144: PetscLogEventBegin(MAT_CholeskyFactorSymbolic,mat,perm,0,0);
3145: (fact->ops->choleskyfactorsymbolic)(fact,mat,perm,info);
3146: PetscLogEventEnd(MAT_CholeskyFactorSymbolic,mat,perm,0,0);
3147: PetscObjectStateIncrease((PetscObject)fact);
3148: return(0);
3149: }
3151: /*@C
3152: MatCholeskyFactorNumeric - Performs numeric Cholesky factorization
3153: of a symmetric matrix. Call this routine after first calling
3154: MatCholeskyFactorSymbolic().
3156: Collective on Mat
3158: Input Parameters:
3159: + fact - the factor matrix obtained with MatGetFactor()
3160: . mat - the initial matrix
3161: . info - options for factorization
3162: - fact - the symbolic factor of mat
3165: Notes:
3166: Most users should employ the simplified KSP interface for linear solvers
3167: instead of working directly with matrix algebra routines such as this.
3168: See, e.g., KSPCreate().
3170: Level: developer
3172: Concepts: matrices^Cholesky numeric factorization
3174: .seealso: MatCholeskyFactorSymbolic(), MatCholeskyFactor(), MatLUFactorNumeric()
3176: Developer Note: fortran interface is not autogenerated as the f90
3177: interface defintion cannot be generated correctly [due to MatFactorInfo]
3179: @*/
3180: PetscErrorCode MatCholeskyFactorNumeric(Mat fact,Mat mat,const MatFactorInfo *info)
3181: {
3189: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
3190: if (!(fact)->ops->choleskyfactornumeric) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s numeric factor Cholesky",((PetscObject)mat)->type_name);
3191: if (mat->rmap->N != (fact)->rmap->N || mat->cmap->N != (fact)->cmap->N) SETERRQ4(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_SIZ,"Mat mat,Mat fact: global dim %D should = %D %D should = %D",mat->rmap->N,(fact)->rmap->N,mat->cmap->N,(fact)->cmap->N);
3192: MatCheckPreallocated(mat,2);
3194: PetscLogEventBegin(MAT_CholeskyFactorNumeric,mat,fact,0,0);
3195: (fact->ops->choleskyfactornumeric)(fact,mat,info);
3196: PetscLogEventEnd(MAT_CholeskyFactorNumeric,mat,fact,0,0);
3197: MatViewFromOptions(fact,NULL,"-mat_factor_view");
3198: PetscObjectStateIncrease((PetscObject)fact);
3199: return(0);
3200: }
3202: /* ----------------------------------------------------------------*/
3203: /*@
3204: MatSolve - Solves A x = b, given a factored matrix.
3206: Neighbor-wise Collective on Mat and Vec
3208: Input Parameters:
3209: + mat - the factored matrix
3210: - b - the right-hand-side vector
3212: Output Parameter:
3213: . x - the result vector
3215: Notes:
3216: The vectors b and x cannot be the same. I.e., one cannot
3217: call MatSolve(A,x,x).
3219: Notes:
3220: Most users should employ the simplified KSP interface for linear solvers
3221: instead of working directly with matrix algebra routines such as this.
3222: See, e.g., KSPCreate().
3224: Level: developer
3226: Concepts: matrices^triangular solves
3228: .seealso: MatSolveAdd(), MatSolveTranspose(), MatSolveTransposeAdd()
3229: @*/
3230: PetscErrorCode MatSolve(Mat mat,Vec b,Vec x)
3231: {
3241: if (x == b) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_IDN,"x and b must be different vectors");
3242: if (mat->cmap->N != x->map->N) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_SIZ,"Mat mat,Vec x: global dim %D %D",mat->cmap->N,x->map->N);
3243: if (mat->rmap->N != b->map->N) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_SIZ,"Mat mat,Vec b: global dim %D %D",mat->rmap->N,b->map->N);
3244: if (mat->rmap->n != b->map->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat mat,Vec b: local dim %D %D",mat->rmap->n,b->map->n);
3245: if (!mat->rmap->N && !mat->cmap->N) return(0);
3246: if (!mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Unfactored matrix");
3247: if (!mat->ops->solve) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
3248: MatCheckPreallocated(mat,1);
3250: PetscLogEventBegin(MAT_Solve,mat,b,x,0);
3251: if (mat->factorerrortype) {
3252: PetscInfo1(mat,"MatFactorError %D\n",mat->factorerrortype);
3253: VecSetInf(x);
3254: } else {
3255: (*mat->ops->solve)(mat,b,x);
3256: }
3257: PetscLogEventEnd(MAT_Solve,mat,b,x,0);
3258: PetscObjectStateIncrease((PetscObject)x);
3259: return(0);
3260: }
3262: static PetscErrorCode MatMatSolve_Basic(Mat A,Mat B,Mat X, PetscBool trans)
3263: {
3265: Vec b,x;
3266: PetscInt m,N,i;
3267: PetscScalar *bb,*xx;
3268: PetscBool flg;
3271: PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,NULL);
3272: if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
3273: PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);
3274: if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");
3276: MatDenseGetArray(B,&bb);
3277: MatDenseGetArray(X,&xx);
3278: MatGetLocalSize(B,&m,NULL); /* number local rows */
3279: MatGetSize(B,NULL,&N); /* total columns in dense matrix */
3280: MatCreateVecs(A,&x,&b);
3281: for (i=0; i<N; i++) {
3282: VecPlaceArray(b,bb + i*m);
3283: VecPlaceArray(x,xx + i*m);
3284: if (trans) {
3285: MatSolveTranspose(A,b,x);
3286: } else {
3287: MatSolve(A,b,x);
3288: }
3289: VecResetArray(x);
3290: VecResetArray(b);
3291: }
3292: VecDestroy(&b);
3293: VecDestroy(&x);
3294: MatDenseRestoreArray(B,&bb);
3295: MatDenseRestoreArray(X,&xx);
3296: return(0);
3297: }
3299: /*@
3300: MatMatSolve - Solves A X = B, given a factored matrix.
3302: Neighbor-wise Collective on Mat
3304: Input Parameters:
3305: + A - the factored matrix
3306: - B - the right-hand-side matrix (dense matrix)
3308: Output Parameter:
3309: . X - the result matrix (dense matrix)
3311: Notes:
3312: The matrices b and x cannot be the same. I.e., one cannot
3313: call MatMatSolve(A,x,x).
3315: Notes:
3316: Most users should usually employ the simplified KSP interface for linear solvers
3317: instead of working directly with matrix algebra routines such as this.
3318: See, e.g., KSPCreate(). However KSP can only solve for one vector (column of X)
3319: at a time.
3321: When using SuperLU_Dist as a parallel solver PETSc will use the SuperLU_Dist functionality to solve multiple right hand sides simultaneously. For MUMPS
3322: it calls a separate solve for each right hand side since MUMPS does not yet support distributed right hand sides.
3324: Since the resulting matrix X must always be dense we do not support sparse representation of the matrix B.
3326: Level: developer
3328: Concepts: matrices^triangular solves
3330: .seealso: MatMatSolveTranspose(), MatLUFactor(), MatCholeskyFactor()
3331: @*/
3332: PetscErrorCode MatMatSolve(Mat A,Mat B,Mat X)
3333: {
3343: if (X == B) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_IDN,"X and B must be different matrices");
3344: if (A->cmap->N != X->rmap->N) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Mat A,Mat X: global dim %D %D",A->cmap->N,X->rmap->N);
3345: if (A->rmap->N != B->rmap->N) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Mat A,Mat B: global dim %D %D",A->rmap->N,B->rmap->N);
3346: if (A->rmap->n != B->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat A,Mat B: local dim %D %D",A->rmap->n,B->rmap->n);
3347: if (X->cmap->N < B->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Solution matrix must have same number of columns as rhs matrix");
3348: if (!A->rmap->N && !A->cmap->N) return(0);
3349: if (!A->factortype) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Unfactored matrix");
3350: MatCheckPreallocated(A,1);
3352: PetscLogEventBegin(MAT_MatSolve,A,B,X,0);
3353: if (!A->ops->matsolve) {
3354: PetscInfo1(A,"Mat type %s using basic MatMatSolve\n",((PetscObject)A)->type_name);
3355: MatMatSolve_Basic(A,B,X,PETSC_FALSE);
3356: } else {
3357: (*A->ops->matsolve)(A,B,X);
3358: }
3359: PetscLogEventEnd(MAT_MatSolve,A,B,X,0);
3360: PetscObjectStateIncrease((PetscObject)X);
3361: return(0);
3362: }
3364: /*@
3365: MatMatSolveTranspose - Solves A^T X = B, given a factored matrix.
3367: Neighbor-wise Collective on Mat
3369: Input Parameters:
3370: + A - the factored matrix
3371: - B - the right-hand-side matrix (dense matrix)
3373: Output Parameter:
3374: . X - the result matrix (dense matrix)
3376: Notes:
3377: The matrices b and x cannot be the same. I.e., one cannot
3378: call MatMatSolveTranspose(A,x,x).
3380: Notes:
3381: Most users should usually employ the simplified KSP interface for linear solvers
3382: instead of working directly with matrix algebra routines such as this.
3383: See, e.g., KSPCreate(). However KSP can only solve for one vector (column of X)
3384: at a time.
3386: When using SuperLU_Dist as a parallel solver PETSc will use the SuperLU_Dist functionality to solve multiple right hand sides simultaneously. For MUMPS
3387: it calls a separate solve for each right hand side since MUMPS does not yet support distributed right hand sides.
3389: Since the resulting matrix X must always be dense we do not support sparse representation of the matrix B.
3391: Level: developer
3393: Concepts: matrices^triangular solves
3395: .seealso: MatMatSolveTranspose(), MatLUFactor(), MatCholeskyFactor()
3396: @*/
3397: PetscErrorCode MatMatSolveTranspose(Mat A,Mat B,Mat X)
3398: {
3408: if (X == B) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_IDN,"X and B must be different matrices");
3409: if (A->cmap->N != X->rmap->N) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Mat A,Mat X: global dim %D %D",A->cmap->N,X->rmap->N);
3410: if (A->rmap->N != B->rmap->N) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Mat A,Mat B: global dim %D %D",A->rmap->N,B->rmap->N);
3411: if (A->rmap->n != B->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat A,Mat B: local dim %D %D",A->rmap->n,B->rmap->n);
3412: if (X->cmap->N < B->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Solution matrix must have same number of columns as rhs matrix");
3413: if (!A->rmap->N && !A->cmap->N) return(0);
3414: if (!A->factortype) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Unfactored matrix");
3415: MatCheckPreallocated(A,1);
3417: PetscLogEventBegin(MAT_MatSolve,A,B,X,0);
3418: if (!A->ops->matsolvetranspose) {
3419: PetscInfo1(A,"Mat type %s using basic MatMatSolveTranspose\n",((PetscObject)A)->type_name);
3420: MatMatSolve_Basic(A,B,X,PETSC_TRUE);
3421: } else {
3422: (*A->ops->matsolvetranspose)(A,B,X);
3423: }
3424: PetscLogEventEnd(MAT_MatSolve,A,B,X,0);
3425: PetscObjectStateIncrease((PetscObject)X);
3426: return(0);
3427: }
3429: /*@
3430: MatForwardSolve - Solves L x = b, given a factored matrix, A = LU, or
3431: U^T*D^(1/2) x = b, given a factored symmetric matrix, A = U^T*D*U,
3433: Neighbor-wise Collective on Mat and Vec
3435: Input Parameters:
3436: + mat - the factored matrix
3437: - b - the right-hand-side vector
3439: Output Parameter:
3440: . x - the result vector
3442: Notes:
3443: MatSolve() should be used for most applications, as it performs
3444: a forward solve followed by a backward solve.
3446: The vectors b and x cannot be the same, i.e., one cannot
3447: call MatForwardSolve(A,x,x).
3449: For matrix in seqsbaij format with block size larger than 1,
3450: the diagonal blocks are not implemented as D = D^(1/2) * D^(1/2) yet.
3451: MatForwardSolve() solves U^T*D y = b, and
3452: MatBackwardSolve() solves U x = y.
3453: Thus they do not provide a symmetric preconditioner.
3455: Most users should employ the simplified KSP interface for linear solvers
3456: instead of working directly with matrix algebra routines such as this.
3457: See, e.g., KSPCreate().
3459: Level: developer
3461: Concepts: matrices^forward solves
3463: .seealso: MatSolve(), MatBackwardSolve()
3464: @*/
3465: PetscErrorCode MatForwardSolve(Mat mat,Vec b,Vec x)
3466: {
3476: if (x == b) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_IDN,"x and b must be different vectors");
3477: if (!mat->ops->forwardsolve) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
3478: if (mat->cmap->N != x->map->N) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_SIZ,"Mat mat,Vec x: global dim %D %D",mat->cmap->N,x->map->N);
3479: if (mat->rmap->N != b->map->N) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_SIZ,"Mat mat,Vec b: global dim %D %D",mat->rmap->N,b->map->N);
3480: if (mat->rmap->n != b->map->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat mat,Vec b: local dim %D %D",mat->rmap->n,b->map->n);
3481: if (!mat->rmap->N && !mat->cmap->N) return(0);
3482: if (!mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Unfactored matrix");
3483: MatCheckPreallocated(mat,1);
3484: PetscLogEventBegin(MAT_ForwardSolve,mat,b,x,0);
3485: (*mat->ops->forwardsolve)(mat,b,x);
3486: PetscLogEventEnd(MAT_ForwardSolve,mat,b,x,0);
3487: PetscObjectStateIncrease((PetscObject)x);
3488: return(0);
3489: }
3491: /*@
3492: MatBackwardSolve - Solves U x = b, given a factored matrix, A = LU.
3493: D^(1/2) U x = b, given a factored symmetric matrix, A = U^T*D*U,
3495: Neighbor-wise Collective on Mat and Vec
3497: Input Parameters:
3498: + mat - the factored matrix
3499: - b - the right-hand-side vector
3501: Output Parameter:
3502: . x - the result vector
3504: Notes:
3505: MatSolve() should be used for most applications, as it performs
3506: a forward solve followed by a backward solve.
3508: The vectors b and x cannot be the same. I.e., one cannot
3509: call MatBackwardSolve(A,x,x).
3511: For matrix in seqsbaij format with block size larger than 1,
3512: the diagonal blocks are not implemented as D = D^(1/2) * D^(1/2) yet.
3513: MatForwardSolve() solves U^T*D y = b, and
3514: MatBackwardSolve() solves U x = y.
3515: Thus they do not provide a symmetric preconditioner.
3517: Most users should employ the simplified KSP interface for linear solvers
3518: instead of working directly with matrix algebra routines such as this.
3519: See, e.g., KSPCreate().
3521: Level: developer
3523: Concepts: matrices^backward solves
3525: .seealso: MatSolve(), MatForwardSolve()
3526: @*/
3527: PetscErrorCode MatBackwardSolve(Mat mat,Vec b,Vec x)
3528: {
3538: if (x == b) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_IDN,"x and b must be different vectors");
3539: if (!mat->ops->backwardsolve) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
3540: if (mat->cmap->N != x->map->N) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_SIZ,"Mat mat,Vec x: global dim %D %D",mat->cmap->N,x->map->N);
3541: if (mat->rmap->N != b->map->N) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_SIZ,"Mat mat,Vec b: global dim %D %D",mat->rmap->N,b->map->N);
3542: if (mat->rmap->n != b->map->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat mat,Vec b: local dim %D %D",mat->rmap->n,b->map->n);
3543: if (!mat->rmap->N && !mat->cmap->N) return(0);
3544: if (!mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Unfactored matrix");
3545: MatCheckPreallocated(mat,1);
3547: PetscLogEventBegin(MAT_BackwardSolve,mat,b,x,0);
3548: (*mat->ops->backwardsolve)(mat,b,x);
3549: PetscLogEventEnd(MAT_BackwardSolve,mat,b,x,0);
3550: PetscObjectStateIncrease((PetscObject)x);
3551: return(0);
3552: }
3554: /*@
3555: MatSolveAdd - Computes x = y + inv(A)*b, given a factored matrix.
3557: Neighbor-wise Collective on Mat and Vec
3559: Input Parameters:
3560: + mat - the factored matrix
3561: . b - the right-hand-side vector
3562: - y - the vector to be added to
3564: Output Parameter:
3565: . x - the result vector
3567: Notes:
3568: The vectors b and x cannot be the same. I.e., one cannot
3569: call MatSolveAdd(A,x,y,x).
3571: Most users should employ the simplified KSP interface for linear solvers
3572: instead of working directly with matrix algebra routines such as this.
3573: See, e.g., KSPCreate().
3575: Level: developer
3577: Concepts: matrices^triangular solves
3579: .seealso: MatSolve(), MatSolveTranspose(), MatSolveTransposeAdd()
3580: @*/
3581: PetscErrorCode MatSolveAdd(Mat mat,Vec b,Vec y,Vec x)
3582: {
3583: PetscScalar one = 1.0;
3584: Vec tmp;
3596: if (x == b) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_IDN,"x and b must be different vectors");
3597: if (mat->cmap->N != x->map->N) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_SIZ,"Mat mat,Vec x: global dim %D %D",mat->cmap->N,x->map->N);
3598: if (mat->rmap->N != b->map->N) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_SIZ,"Mat mat,Vec b: global dim %D %D",mat->rmap->N,b->map->N);
3599: if (mat->rmap->N != y->map->N) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_SIZ,"Mat mat,Vec y: global dim %D %D",mat->rmap->N,y->map->N);
3600: if (mat->rmap->n != b->map->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat mat,Vec b: local dim %D %D",mat->rmap->n,b->map->n);
3601: if (x->map->n != y->map->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Vec x,Vec y: local dim %D %D",x->map->n,y->map->n);
3602: if (!mat->rmap->N && !mat->cmap->N) return(0);
3603: if (!mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Unfactored matrix");
3604: MatCheckPreallocated(mat,1);
3606: PetscLogEventBegin(MAT_SolveAdd,mat,b,x,y);
3607: if (mat->ops->solveadd) {
3608: (*mat->ops->solveadd)(mat,b,y,x);
3609: } else {
3610: /* do the solve then the add manually */
3611: if (x != y) {
3612: MatSolve(mat,b,x);
3613: VecAXPY(x,one,y);
3614: } else {
3615: VecDuplicate(x,&tmp);
3616: PetscLogObjectParent((PetscObject)mat,(PetscObject)tmp);
3617: VecCopy(x,tmp);
3618: MatSolve(mat,b,x);
3619: VecAXPY(x,one,tmp);
3620: VecDestroy(&tmp);
3621: }
3622: }
3623: PetscLogEventEnd(MAT_SolveAdd,mat,b,x,y);
3624: PetscObjectStateIncrease((PetscObject)x);
3625: return(0);
3626: }
3628: /*@
3629: MatSolveTranspose - Solves A' x = b, given a factored matrix.
3631: Neighbor-wise Collective on Mat and Vec
3633: Input Parameters:
3634: + mat - the factored matrix
3635: - b - the right-hand-side vector
3637: Output Parameter:
3638: . x - the result vector
3640: Notes:
3641: The vectors b and x cannot be the same. I.e., one cannot
3642: call MatSolveTranspose(A,x,x).
3644: Most users should employ the simplified KSP interface for linear solvers
3645: instead of working directly with matrix algebra routines such as this.
3646: See, e.g., KSPCreate().
3648: Level: developer
3650: Concepts: matrices^triangular solves
3652: .seealso: MatSolve(), MatSolveAdd(), MatSolveTransposeAdd()
3653: @*/
3654: PetscErrorCode MatSolveTranspose(Mat mat,Vec b,Vec x)
3655: {
3665: if (x == b) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_IDN,"x and b must be different vectors");
3666: if (!mat->ops->solvetranspose) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Matrix type %s",((PetscObject)mat)->type_name);
3667: if (mat->rmap->N != x->map->N) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_SIZ,"Mat mat,Vec x: global dim %D %D",mat->rmap->N,x->map->N);
3668: if (mat->cmap->N != b->map->N) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_SIZ,"Mat mat,Vec b: global dim %D %D",mat->cmap->N,b->map->N);
3669: if (!mat->rmap->N && !mat->cmap->N) return(0);
3670: if (!mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Unfactored matrix");
3671: MatCheckPreallocated(mat,1);
3672: PetscLogEventBegin(MAT_SolveTranspose,mat,b,x,0);
3673: if (mat->factorerrortype) {
3674: PetscInfo1(mat,"MatFactorError %D\n",mat->factorerrortype);
3675: VecSetInf(x);
3676: } else {
3677: (*mat->ops->solvetranspose)(mat,b,x);
3678: }
3679: PetscLogEventEnd(MAT_SolveTranspose,mat,b,x,0);
3680: PetscObjectStateIncrease((PetscObject)x);
3681: return(0);
3682: }
3684: /*@
3685: MatSolveTransposeAdd - Computes x = y + inv(Transpose(A)) b, given a
3686: factored matrix.
3688: Neighbor-wise Collective on Mat and Vec
3690: Input Parameters:
3691: + mat - the factored matrix
3692: . b - the right-hand-side vector
3693: - y - the vector to be added to
3695: Output Parameter:
3696: . x - the result vector
3698: Notes:
3699: The vectors b and x cannot be the same. I.e., one cannot
3700: call MatSolveTransposeAdd(A,x,y,x).
3702: Most users should employ the simplified KSP interface for linear solvers
3703: instead of working directly with matrix algebra routines such as this.
3704: See, e.g., KSPCreate().
3706: Level: developer
3708: Concepts: matrices^triangular solves
3710: .seealso: MatSolve(), MatSolveAdd(), MatSolveTranspose()
3711: @*/
3712: PetscErrorCode MatSolveTransposeAdd(Mat mat,Vec b,Vec y,Vec x)
3713: {
3714: PetscScalar one = 1.0;
3716: Vec tmp;
3727: if (x == b) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_IDN,"x and b must be different vectors");
3728: if (mat->rmap->N != x->map->N) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_SIZ,"Mat mat,Vec x: global dim %D %D",mat->rmap->N,x->map->N);
3729: if (mat->cmap->N != b->map->N) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_SIZ,"Mat mat,Vec b: global dim %D %D",mat->cmap->N,b->map->N);
3730: if (mat->cmap->N != y->map->N) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_SIZ,"Mat mat,Vec y: global dim %D %D",mat->cmap->N,y->map->N);
3731: if (x->map->n != y->map->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Vec x,Vec y: local dim %D %D",x->map->n,y->map->n);
3732: if (!mat->rmap->N && !mat->cmap->N) return(0);
3733: if (!mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Unfactored matrix");
3734: MatCheckPreallocated(mat,1);
3736: PetscLogEventBegin(MAT_SolveTransposeAdd,mat,b,x,y);
3737: if (mat->ops->solvetransposeadd) {
3738: if (mat->factorerrortype) {
3739: PetscInfo1(mat,"MatFactorError %D\n",mat->factorerrortype);
3740: VecSetInf(x);
3741: } else {
3742: (*mat->ops->solvetransposeadd)(mat,b,y,x);
3743: }
3744: } else {
3745: /* do the solve then the add manually */
3746: if (x != y) {
3747: MatSolveTranspose(mat,b,x);
3748: VecAXPY(x,one,y);
3749: } else {
3750: VecDuplicate(x,&tmp);
3751: PetscLogObjectParent((PetscObject)mat,(PetscObject)tmp);
3752: VecCopy(x,tmp);
3753: MatSolveTranspose(mat,b,x);
3754: VecAXPY(x,one,tmp);
3755: VecDestroy(&tmp);
3756: }
3757: }
3758: PetscLogEventEnd(MAT_SolveTransposeAdd,mat,b,x,y);
3759: PetscObjectStateIncrease((PetscObject)x);
3760: return(0);
3761: }
3762: /* ----------------------------------------------------------------*/
3764: /*@
3765: MatSOR - Computes relaxation (SOR, Gauss-Seidel) sweeps.
3767: Neighbor-wise Collective on Mat and Vec
3769: Input Parameters:
3770: + mat - the matrix
3771: . b - the right hand side
3772: . omega - the relaxation factor
3773: . flag - flag indicating the type of SOR (see below)
3774: . shift - diagonal shift
3775: . its - the number of iterations
3776: - lits - the number of local iterations
3778: Output Parameters:
3779: . x - the solution (can contain an initial guess, use option SOR_ZERO_INITIAL_GUESS to indicate no guess)
3781: SOR Flags:
3782: . SOR_FORWARD_SWEEP - forward SOR
3783: . SOR_BACKWARD_SWEEP - backward SOR
3784: . SOR_SYMMETRIC_SWEEP - SSOR (symmetric SOR)
3785: . SOR_LOCAL_FORWARD_SWEEP - local forward SOR
3786: . SOR_LOCAL_BACKWARD_SWEEP - local forward SOR
3787: . SOR_LOCAL_SYMMETRIC_SWEEP - local SSOR
3788: . SOR_APPLY_UPPER, SOR_APPLY_LOWER - applies
3789: upper/lower triangular part of matrix to
3790: vector (with omega)
3791: . SOR_ZERO_INITIAL_GUESS - zero initial guess
3793: Notes:
3794: SOR_LOCAL_FORWARD_SWEEP, SOR_LOCAL_BACKWARD_SWEEP, and
3795: SOR_LOCAL_SYMMETRIC_SWEEP perform separate independent smoothings
3796: on each processor.
3798: Application programmers will not generally use MatSOR() directly,
3799: but instead will employ the KSP/PC interface.
3801: Notes: for BAIJ, SBAIJ, and AIJ matrices with Inodes this does a block SOR smoothing, otherwise it does a pointwise smoothing
3803: Notes for Advanced Users:
3804: The flags are implemented as bitwise inclusive or operations.
3805: For example, use (SOR_ZERO_INITIAL_GUESS | SOR_SYMMETRIC_SWEEP)
3806: to specify a zero initial guess for SSOR.
3808: Most users should employ the simplified KSP interface for linear solvers
3809: instead of working directly with matrix algebra routines such as this.
3810: See, e.g., KSPCreate().
3812: Vectors x and b CANNOT be the same
3814: Developer Note: We should add block SOR support for AIJ matrices with block size set to great than one and no inodes
3816: Level: developer
3818: Concepts: matrices^relaxation
3819: Concepts: matrices^SOR
3820: Concepts: matrices^Gauss-Seidel
3822: @*/
3823: PetscErrorCode MatSOR(Mat mat,Vec b,PetscReal omega,MatSORType flag,PetscReal shift,PetscInt its,PetscInt lits,Vec x)
3824: {
3834: if (!mat->ops->sor) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
3835: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
3836: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
3837: if (mat->cmap->N != x->map->N) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_SIZ,"Mat mat,Vec x: global dim %D %D",mat->cmap->N,x->map->N);
3838: if (mat->rmap->N != b->map->N) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_SIZ,"Mat mat,Vec b: global dim %D %D",mat->rmap->N,b->map->N);
3839: if (mat->rmap->n != b->map->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat mat,Vec b: local dim %D %D",mat->rmap->n,b->map->n);
3840: if (its <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D positive",its);
3841: if (lits <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires local its %D positive",lits);
3842: if (b == x) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_IDN,"b and x vector cannot be the same");
3844: MatCheckPreallocated(mat,1);
3845: PetscLogEventBegin(MAT_SOR,mat,b,x,0);
3846: ierr =(*mat->ops->sor)(mat,b,omega,flag,shift,its,lits,x);
3847: PetscLogEventEnd(MAT_SOR,mat,b,x,0);
3848: PetscObjectStateIncrease((PetscObject)x);
3849: return(0);
3850: }
3852: /*
3853: Default matrix copy routine.
3854: */
3855: PetscErrorCode MatCopy_Basic(Mat A,Mat B,MatStructure str)
3856: {
3857: PetscErrorCode ierr;
3858: PetscInt i,rstart = 0,rend = 0,nz;
3859: const PetscInt *cwork;
3860: const PetscScalar *vwork;
3863: if (B->assembled) {
3864: MatZeroEntries(B);
3865: }
3866: MatGetOwnershipRange(A,&rstart,&rend);
3867: for (i=rstart; i<rend; i++) {
3868: MatGetRow(A,i,&nz,&cwork,&vwork);
3869: MatSetValues(B,1,&i,nz,cwork,vwork,INSERT_VALUES);
3870: MatRestoreRow(A,i,&nz,&cwork,&vwork);
3871: }
3872: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
3873: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
3874: return(0);
3875: }
3877: /*@
3878: MatCopy - Copys a matrix to another matrix.
3880: Collective on Mat
3882: Input Parameters:
3883: + A - the matrix
3884: - str - SAME_NONZERO_PATTERN or DIFFERENT_NONZERO_PATTERN
3886: Output Parameter:
3887: . B - where the copy is put
3889: Notes:
3890: If you use SAME_NONZERO_PATTERN then the two matrices had better have the
3891: same nonzero pattern or the routine will crash.
3893: MatCopy() copies the matrix entries of a matrix to another existing
3894: matrix (after first zeroing the second matrix). A related routine is
3895: MatConvert(), which first creates a new matrix and then copies the data.
3897: Level: intermediate
3899: Concepts: matrices^copying
3901: .seealso: MatConvert(), MatDuplicate()
3903: @*/
3904: PetscErrorCode MatCopy(Mat A,Mat B,MatStructure str)
3905: {
3907: PetscInt i;
3915: MatCheckPreallocated(B,2);
3916: if (!A->assembled) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
3917: if (A->factortype) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
3918: if (A->rmap->N != B->rmap->N || A->cmap->N != B->cmap->N) SETERRQ4(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Mat A,Mat B: global dim (%D,%D) (%D,%D)",A->rmap->N,B->rmap->N,A->cmap->N,B->cmap->N);
3919: MatCheckPreallocated(A,1);
3920: if (A == B) return(0);
3922: PetscLogEventBegin(MAT_Copy,A,B,0,0);
3923: if (A->ops->copy) {
3924: (*A->ops->copy)(A,B,str);
3925: } else { /* generic conversion */
3926: MatCopy_Basic(A,B,str);
3927: }
3929: B->stencil.dim = A->stencil.dim;
3930: B->stencil.noc = A->stencil.noc;
3931: for (i=0; i<=A->stencil.dim; i++) {
3932: B->stencil.dims[i] = A->stencil.dims[i];
3933: B->stencil.starts[i] = A->stencil.starts[i];
3934: }
3936: PetscLogEventEnd(MAT_Copy,A,B,0,0);
3937: PetscObjectStateIncrease((PetscObject)B);
3938: return(0);
3939: }
3941: /*@C
3942: MatConvert - Converts a matrix to another matrix, either of the same
3943: or different type.
3945: Collective on Mat
3947: Input Parameters:
3948: + mat - the matrix
3949: . newtype - new matrix type. Use MATSAME to create a new matrix of the
3950: same type as the original matrix.
3951: - reuse - denotes if the destination matrix is to be created or reused.
3952: Use MAT_INPLACE_MATRIX for inplace conversion (that is when you want the input mat to be changed to contain the matrix in the new format), otherwise use
3953: MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX (can only be used after the first call was made with MAT_INITIAL_MATRIX, causes the matrix space in M to be reused).
3955: Output Parameter:
3956: . M - pointer to place new matrix
3958: Notes:
3959: MatConvert() first creates a new matrix and then copies the data from
3960: the first matrix. A related routine is MatCopy(), which copies the matrix
3961: entries of one matrix to another already existing matrix context.
3963: Cannot be used to convert a sequential matrix to parallel or parallel to sequential,
3964: the MPI communicator of the generated matrix is always the same as the communicator
3965: of the input matrix.
3967: Level: intermediate
3969: Concepts: matrices^converting between storage formats
3971: .seealso: MatCopy(), MatDuplicate()
3972: @*/
3973: PetscErrorCode MatConvert(Mat mat, MatType newtype,MatReuse reuse,Mat *M)
3974: {
3976: PetscBool sametype,issame,flg;
3977: char convname[256],mtype[256];
3978: Mat B;
3984: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
3985: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
3986: MatCheckPreallocated(mat,1);
3988: PetscOptionsGetString(((PetscObject)mat)->options,((PetscObject)mat)->prefix,"-matconvert_type",mtype,256,&flg);
3989: if (flg) {
3990: newtype = mtype;
3991: }
3992: PetscObjectTypeCompare((PetscObject)mat,newtype,&sametype);
3993: PetscStrcmp(newtype,"same",&issame);
3994: if ((reuse == MAT_INPLACE_MATRIX) && (mat != *M)) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"MAT_INPLACE_MATRIX requires same input and output matrix");
3995: if ((reuse == MAT_REUSE_MATRIX) && (mat == *M)) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"MAT_REUSE_MATRIX means reuse matrix in final argument, perhaps you mean MAT_INPLACE_MATRIX");
3997: if ((reuse == MAT_INPLACE_MATRIX) && (issame || sametype)) return(0);
3999: if ((sametype || issame) && (reuse==MAT_INITIAL_MATRIX) && mat->ops->duplicate) {
4000: (*mat->ops->duplicate)(mat,MAT_COPY_VALUES,M);
4001: } else {
4002: PetscErrorCode (*conv)(Mat, MatType,MatReuse,Mat*)=NULL;
4003: const char *prefix[3] = {"seq","mpi",""};
4004: PetscInt i;
4005: /*
4006: Order of precedence:
4007: 1) See if a specialized converter is known to the current matrix.
4008: 2) See if a specialized converter is known to the desired matrix class.
4009: 3) See if a good general converter is registered for the desired class
4010: (as of 6/27/03 only MATMPIADJ falls into this category).
4011: 4) See if a good general converter is known for the current matrix.
4012: 5) Use a really basic converter.
4013: */
4015: /* 1) See if a specialized converter is known to the current matrix and the desired class */
4016: for (i=0; i<3; i++) {
4017: PetscStrncpy(convname,"MatConvert_",sizeof(convname));
4018: PetscStrlcat(convname,((PetscObject)mat)->type_name,sizeof(convname));
4019: PetscStrlcat(convname,"_",sizeof(convname));
4020: PetscStrlcat(convname,prefix[i],sizeof(convname));
4021: PetscStrlcat(convname,issame ? ((PetscObject)mat)->type_name : newtype,sizeof(convname));
4022: PetscStrlcat(convname,"_C",sizeof(convname));
4023: PetscObjectQueryFunction((PetscObject)mat,convname,&conv);
4024: if (conv) goto foundconv;
4025: }
4027: /* 2) See if a specialized converter is known to the desired matrix class. */
4028: MatCreate(PetscObjectComm((PetscObject)mat),&B);
4029: MatSetSizes(B,mat->rmap->n,mat->cmap->n,mat->rmap->N,mat->cmap->N);
4030: MatSetType(B,newtype);
4031: for (i=0; i<3; i++) {
4032: PetscStrncpy(convname,"MatConvert_",sizeof(convname));
4033: PetscStrlcat(convname,((PetscObject)mat)->type_name,sizeof(convname));
4034: PetscStrlcat(convname,"_",sizeof(convname));
4035: PetscStrlcat(convname,prefix[i],sizeof(convname));
4036: PetscStrlcat(convname,newtype,sizeof(convname));
4037: PetscStrlcat(convname,"_C",sizeof(convname));
4038: PetscObjectQueryFunction((PetscObject)B,convname,&conv);
4039: if (conv) {
4040: MatDestroy(&B);
4041: goto foundconv;
4042: }
4043: }
4045: /* 3) See if a good general converter is registered for the desired class */
4046: conv = B->ops->convertfrom;
4047: MatDestroy(&B);
4048: if (conv) goto foundconv;
4050: /* 4) See if a good general converter is known for the current matrix */
4051: if (mat->ops->convert) {
4052: conv = mat->ops->convert;
4053: }
4054: if (conv) goto foundconv;
4056: /* 5) Use a really basic converter. */
4057: conv = MatConvert_Basic;
4059: foundconv:
4060: PetscLogEventBegin(MAT_Convert,mat,0,0,0);
4061: (*conv)(mat,newtype,reuse,M);
4062: if (mat->rmap->mapping && mat->cmap->mapping && !(*M)->rmap->mapping && !(*M)->cmap->mapping) {
4063: /* the block sizes must be same if the mappings are copied over */
4064: (*M)->rmap->bs = mat->rmap->bs;
4065: (*M)->cmap->bs = mat->cmap->bs;
4066: PetscObjectReference((PetscObject)mat->rmap->mapping);
4067: PetscObjectReference((PetscObject)mat->cmap->mapping);
4068: (*M)->rmap->mapping = mat->rmap->mapping;
4069: (*M)->cmap->mapping = mat->cmap->mapping;
4070: }
4071: (*M)->stencil.dim = mat->stencil.dim;
4072: (*M)->stencil.noc = mat->stencil.noc;
4073: for (i=0; i<=mat->stencil.dim; i++) {
4074: (*M)->stencil.dims[i] = mat->stencil.dims[i];
4075: (*M)->stencil.starts[i] = mat->stencil.starts[i];
4076: }
4077: PetscLogEventEnd(MAT_Convert,mat,0,0,0);
4078: }
4079: PetscObjectStateIncrease((PetscObject)*M);
4081: /* Copy Mat options */
4082: if (mat->symmetric) {MatSetOption(*M,MAT_SYMMETRIC,PETSC_TRUE);}
4083: if (mat->hermitian) {MatSetOption(*M,MAT_HERMITIAN,PETSC_TRUE);}
4084: return(0);
4085: }
4087: /*@C
4088: MatFactorGetSolverType - Returns name of the package providing the factorization routines
4090: Not Collective
4092: Input Parameter:
4093: . mat - the matrix, must be a factored matrix
4095: Output Parameter:
4096: . type - the string name of the package (do not free this string)
4098: Notes:
4099: In Fortran you pass in a empty string and the package name will be copied into it.
4100: (Make sure the string is long enough)
4102: Level: intermediate
4104: .seealso: MatCopy(), MatDuplicate(), MatGetFactorAvailable(), MatGetFactor()
4105: @*/
4106: PetscErrorCode MatFactorGetSolverType(Mat mat, MatSolverType *type)
4107: {
4108: PetscErrorCode ierr, (*conv)(Mat,MatSolverType*);
4113: if (!mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
4114: PetscObjectQueryFunction((PetscObject)mat,"MatFactorGetSolverType_C",&conv);
4115: if (!conv) {
4116: *type = MATSOLVERPETSC;
4117: } else {
4118: (*conv)(mat,type);
4119: }
4120: return(0);
4121: }
4123: typedef struct _MatSolverTypeForSpecifcType* MatSolverTypeForSpecifcType;
4124: struct _MatSolverTypeForSpecifcType {
4125: MatType mtype;
4126: PetscErrorCode (*getfactor[4])(Mat,MatFactorType,Mat*);
4127: MatSolverTypeForSpecifcType next;
4128: };
4130: typedef struct _MatSolverTypeHolder* MatSolverTypeHolder;
4131: struct _MatSolverTypeHolder {
4132: char *name;
4133: MatSolverTypeForSpecifcType handlers;
4134: MatSolverTypeHolder next;
4135: };
4137: static MatSolverTypeHolder MatSolverTypeHolders = NULL;
4139: /*@C
4140: MatSolvePackageRegister - Registers a MatSolverType that works for a particular matrix type
4142: Input Parameters:
4143: + package - name of the package, for example petsc or superlu
4144: . mtype - the matrix type that works with this package
4145: . ftype - the type of factorization supported by the package
4146: - getfactor - routine that will create the factored matrix ready to be used
4148: Level: intermediate
4150: .seealso: MatCopy(), MatDuplicate(), MatGetFactorAvailable()
4151: @*/
4152: PetscErrorCode MatSolverTypeRegister(MatSolverType package,const MatType mtype,MatFactorType ftype,PetscErrorCode (*getfactor)(Mat,MatFactorType,Mat*))
4153: {
4154: PetscErrorCode ierr;
4155: MatSolverTypeHolder next = MatSolverTypeHolders,prev;
4156: PetscBool flg;
4157: MatSolverTypeForSpecifcType inext,iprev = NULL;
4160: if (!next) {
4161: PetscNew(&MatSolverTypeHolders);
4162: PetscStrallocpy(package,&MatSolverTypeHolders->name);
4163: PetscNew(&MatSolverTypeHolders->handlers);
4164: PetscStrallocpy(mtype,(char **)&MatSolverTypeHolders->handlers->mtype);
4165: MatSolverTypeHolders->handlers->getfactor[(int)ftype-1] = getfactor;
4166: return(0);
4167: }
4168: while (next) {
4169: PetscStrcasecmp(package,next->name,&flg);
4170: if (flg) {
4171: if (!next->handlers) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatSolverTypeHolder is missing handlers");
4172: inext = next->handlers;
4173: while (inext) {
4174: PetscStrcasecmp(mtype,inext->mtype,&flg);
4175: if (flg) {
4176: inext->getfactor[(int)ftype-1] = getfactor;
4177: return(0);
4178: }
4179: iprev = inext;
4180: inext = inext->next;
4181: }
4182: PetscNew(&iprev->next);
4183: PetscStrallocpy(mtype,(char **)&iprev->next->mtype);
4184: iprev->next->getfactor[(int)ftype-1] = getfactor;
4185: return(0);
4186: }
4187: prev = next;
4188: next = next->next;
4189: }
4190: PetscNew(&prev->next);
4191: PetscStrallocpy(package,&prev->next->name);
4192: PetscNew(&prev->next->handlers);
4193: PetscStrallocpy(mtype,(char **)&prev->next->handlers->mtype);
4194: prev->next->handlers->getfactor[(int)ftype-1] = getfactor;
4195: return(0);
4196: }
4198: /*@C
4199: MatSolvePackageGet - Get's the function that creates the factor matrix if it exist
4201: Input Parameters:
4202: + package - name of the package, for example petsc or superlu
4203: . ftype - the type of factorization supported by the package
4204: - mtype - the matrix type that works with this package
4206: Output Parameters:
4207: + foundpackage - PETSC_TRUE if the package was registered
4208: . foundmtype - PETSC_TRUE if the package supports the requested mtype
4209: - getfactor - routine that will create the factored matrix ready to be used or NULL if not found
4211: Level: intermediate
4213: .seealso: MatCopy(), MatDuplicate(), MatGetFactorAvailable()
4214: @*/
4215: PetscErrorCode MatSolverTypeGet(MatSolverType package,const MatType mtype,MatFactorType ftype,PetscBool *foundpackage,PetscBool *foundmtype,PetscErrorCode (**getfactor)(Mat,MatFactorType,Mat*))
4216: {
4217: PetscErrorCode ierr;
4218: MatSolverTypeHolder next = MatSolverTypeHolders;
4219: PetscBool flg;
4220: MatSolverTypeForSpecifcType inext;
4223: if (foundpackage) *foundpackage = PETSC_FALSE;
4224: if (foundmtype) *foundmtype = PETSC_FALSE;
4225: if (getfactor) *getfactor = NULL;
4227: if (package) {
4228: while (next) {
4229: PetscStrcasecmp(package,next->name,&flg);
4230: if (flg) {
4231: if (foundpackage) *foundpackage = PETSC_TRUE;
4232: inext = next->handlers;
4233: while (inext) {
4234: PetscStrbeginswith(mtype,inext->mtype,&flg);
4235: if (flg) {
4236: if (foundmtype) *foundmtype = PETSC_TRUE;
4237: if (getfactor) *getfactor = inext->getfactor[(int)ftype-1];
4238: return(0);
4239: }
4240: inext = inext->next;
4241: }
4242: }
4243: next = next->next;
4244: }
4245: } else {
4246: while (next) {
4247: inext = next->handlers;
4248: while (inext) {
4249: PetscStrbeginswith(mtype,inext->mtype,&flg);
4250: if (flg && inext->getfactor[(int)ftype-1]) {
4251: if (foundpackage) *foundpackage = PETSC_TRUE;
4252: if (foundmtype) *foundmtype = PETSC_TRUE;
4253: if (getfactor) *getfactor = inext->getfactor[(int)ftype-1];
4254: return(0);
4255: }
4256: inext = inext->next;
4257: }
4258: next = next->next;
4259: }
4260: }
4261: return(0);
4262: }
4264: PetscErrorCode MatSolverTypeDestroy(void)
4265: {
4266: PetscErrorCode ierr;
4267: MatSolverTypeHolder next = MatSolverTypeHolders,prev;
4268: MatSolverTypeForSpecifcType inext,iprev;
4271: while (next) {
4272: PetscFree(next->name);
4273: inext = next->handlers;
4274: while (inext) {
4275: PetscFree(inext->mtype);
4276: iprev = inext;
4277: inext = inext->next;
4278: PetscFree(iprev);
4279: }
4280: prev = next;
4281: next = next->next;
4282: PetscFree(prev);
4283: }
4284: MatSolverTypeHolders = NULL;
4285: return(0);
4286: }
4288: /*@C
4289: MatGetFactor - Returns a matrix suitable to calls to MatXXFactorSymbolic()
4291: Collective on Mat
4293: Input Parameters:
4294: + mat - the matrix
4295: . type - name of solver type, for example, superlu, petsc (to use PETSc's default)
4296: - ftype - factor type, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ICC, MAT_FACTOR_ILU,
4298: Output Parameters:
4299: . f - the factor matrix used with MatXXFactorSymbolic() calls
4301: Notes:
4302: Some PETSc matrix formats have alternative solvers available that are contained in alternative packages
4303: such as pastix, superlu, mumps etc.
4305: PETSc must have been ./configure to use the external solver, using the option --download-package
4307: Level: intermediate
4309: .seealso: MatCopy(), MatDuplicate(), MatGetFactorAvailable()
4310: @*/
4311: PetscErrorCode MatGetFactor(Mat mat, MatSolverType type,MatFactorType ftype,Mat *f)
4312: {
4313: PetscErrorCode ierr,(*conv)(Mat,MatFactorType,Mat*);
4314: PetscBool foundpackage,foundmtype;
4320: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
4321: MatCheckPreallocated(mat,1);
4323: MatSolverTypeGet(type,((PetscObject)mat)->type_name,ftype,&foundpackage,&foundmtype,&conv);
4324: if (!foundpackage) {
4325: if (type) {
4326: SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_MISSING_FACTOR,"Could not locate solver package %s. Perhaps you must ./configure with --download-%s",type,type);
4327: } else {
4328: SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_MISSING_FACTOR,"Could not locate a solver package. Perhaps you must ./configure with --download-<package>");
4329: }
4330: }
4332: if (!foundmtype) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_MISSING_FACTOR,"MatSolverType %s does not support matrix type %s",type,((PetscObject)mat)->type_name);
4333: if (!conv) SETERRQ3(PetscObjectComm((PetscObject)mat),PETSC_ERR_MISSING_FACTOR,"MatSolverType %s does not support factorization type %s for matrix type %s",type,MatFactorTypes[ftype],((PetscObject)mat)->type_name);
4335: #if defined(PETSC_USE_COMPLEX)
4336: if (mat->hermitian && !mat->symmetric && (ftype == MAT_FACTOR_CHOLESKY||ftype == MAT_FACTOR_ICC)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Hermitian CHOLESKY or ICC Factor is not supported");
4337: #endif
4339: (*conv)(mat,ftype,f);
4340: return(0);
4341: }
4343: /*@C
4344: MatGetFactorAvailable - Returns a a flag if matrix supports particular package and factor type
4346: Not Collective
4348: Input Parameters:
4349: + mat - the matrix
4350: . type - name of solver type, for example, superlu, petsc (to use PETSc's default)
4351: - ftype - factor type, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ICC, MAT_FACTOR_ILU,
4353: Output Parameter:
4354: . flg - PETSC_TRUE if the factorization is available
4356: Notes:
4357: Some PETSc matrix formats have alternative solvers available that are contained in alternative packages
4358: such as pastix, superlu, mumps etc.
4360: PETSc must have been ./configure to use the external solver, using the option --download-package
4362: Level: intermediate
4364: .seealso: MatCopy(), MatDuplicate(), MatGetFactor()
4365: @*/
4366: PetscErrorCode MatGetFactorAvailable(Mat mat, MatSolverType type,MatFactorType ftype,PetscBool *flg)
4367: {
4368: PetscErrorCode ierr, (*gconv)(Mat,MatFactorType,Mat*);
4374: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
4375: MatCheckPreallocated(mat,1);
4377: *flg = PETSC_FALSE;
4378: MatSolverTypeGet(type,((PetscObject)mat)->type_name,ftype,NULL,NULL,&gconv);
4379: if (gconv) {
4380: *flg = PETSC_TRUE;
4381: }
4382: return(0);
4383: }
4385: #include <petscdmtypes.h>
4387: /*@
4388: MatDuplicate - Duplicates a matrix including the non-zero structure.
4390: Collective on Mat
4392: Input Parameters:
4393: + mat - the matrix
4394: - op - One of MAT_DO_NOT_COPY_VALUES, MAT_COPY_VALUES, or MAT_SHARE_NONZERO_PATTERN.
4395: See the manual page for MatDuplicateOption for an explanation of these options.
4397: Output Parameter:
4398: . M - pointer to place new matrix
4400: Level: intermediate
4402: Concepts: matrices^duplicating
4404: Notes: You cannot change the nonzero pattern for the parent or child matrix if you use MAT_SHARE_NONZERO_PATTERN.
4406: .seealso: MatCopy(), MatConvert(), MatDuplicateOption
4407: @*/
4408: PetscErrorCode MatDuplicate(Mat mat,MatDuplicateOption op,Mat *M)
4409: {
4411: Mat B;
4412: PetscInt i;
4413: DM dm;
4414: void (*viewf)(void);
4420: if (op == MAT_COPY_VALUES && !mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"MAT_COPY_VALUES not allowed for unassembled matrix");
4421: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
4422: MatCheckPreallocated(mat,1);
4424: *M = 0;
4425: if (!mat->ops->duplicate) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Not written for this matrix type");
4426: PetscLogEventBegin(MAT_Convert,mat,0,0,0);
4427: (*mat->ops->duplicate)(mat,op,M);
4428: B = *M;
4430: MatGetOperation(mat,MATOP_VIEW,&viewf);
4431: if (viewf) {
4432: MatSetOperation(B,MATOP_VIEW,viewf);
4433: }
4435: B->stencil.dim = mat->stencil.dim;
4436: B->stencil.noc = mat->stencil.noc;
4437: for (i=0; i<=mat->stencil.dim; i++) {
4438: B->stencil.dims[i] = mat->stencil.dims[i];
4439: B->stencil.starts[i] = mat->stencil.starts[i];
4440: }
4442: B->nooffproczerorows = mat->nooffproczerorows;
4443: B->nooffprocentries = mat->nooffprocentries;
4445: PetscObjectQuery((PetscObject) mat, "__PETSc_dm", (PetscObject*) &dm);
4446: if (dm) {
4447: PetscObjectCompose((PetscObject) B, "__PETSc_dm", (PetscObject) dm);
4448: }
4449: PetscLogEventEnd(MAT_Convert,mat,0,0,0);
4450: PetscObjectStateIncrease((PetscObject)B);
4451: return(0);
4452: }
4454: /*@
4455: MatGetDiagonal - Gets the diagonal of a matrix.
4457: Logically Collective on Mat and Vec
4459: Input Parameters:
4460: + mat - the matrix
4461: - v - the vector for storing the diagonal
4463: Output Parameter:
4464: . v - the diagonal of the matrix
4466: Level: intermediate
4468: Note:
4469: Currently only correct in parallel for square matrices.
4471: Concepts: matrices^accessing diagonals
4473: .seealso: MatGetRow(), MatCreateSubMatrices(), MatCreateSubMatrix(), MatGetRowMaxAbs()
4474: @*/
4475: PetscErrorCode MatGetDiagonal(Mat mat,Vec v)
4476: {
4483: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
4484: if (!mat->ops->getdiagonal) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
4485: MatCheckPreallocated(mat,1);
4487: (*mat->ops->getdiagonal)(mat,v);
4488: PetscObjectStateIncrease((PetscObject)v);
4489: return(0);
4490: }
4492: /*@C
4493: MatGetRowMin - Gets the minimum value (of the real part) of each
4494: row of the matrix
4496: Logically Collective on Mat and Vec
4498: Input Parameters:
4499: . mat - the matrix
4501: Output Parameter:
4502: + v - the vector for storing the maximums
4503: - idx - the indices of the column found for each row (optional)
4505: Level: intermediate
4507: Notes: The result of this call are the same as if one converted the matrix to dense format
4508: and found the minimum value in each row (i.e. the implicit zeros are counted as zeros).
4510: This code is only implemented for a couple of matrix formats.
4512: Concepts: matrices^getting row maximums
4514: .seealso: MatGetDiagonal(), MatCreateSubMatrices(), MatCreateSubMatrix(), MatGetRowMaxAbs(),
4515: MatGetRowMax()
4516: @*/
4517: PetscErrorCode MatGetRowMin(Mat mat,Vec v,PetscInt idx[])
4518: {
4525: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
4526: if (!mat->ops->getrowmax) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
4527: MatCheckPreallocated(mat,1);
4529: (*mat->ops->getrowmin)(mat,v,idx);
4530: PetscObjectStateIncrease((PetscObject)v);
4531: return(0);
4532: }
4534: /*@C
4535: MatGetRowMinAbs - Gets the minimum value (in absolute value) of each
4536: row of the matrix
4538: Logically Collective on Mat and Vec
4540: Input Parameters:
4541: . mat - the matrix
4543: Output Parameter:
4544: + v - the vector for storing the minimums
4545: - idx - the indices of the column found for each row (or NULL if not needed)
4547: Level: intermediate
4549: Notes: if a row is completely empty or has only 0.0 values then the idx[] value for that
4550: row is 0 (the first column).
4552: This code is only implemented for a couple of matrix formats.
4554: Concepts: matrices^getting row maximums
4556: .seealso: MatGetDiagonal(), MatCreateSubMatrices(), MatCreateSubMatrix(), MatGetRowMax(), MatGetRowMaxAbs(), MatGetRowMin()
4557: @*/
4558: PetscErrorCode MatGetRowMinAbs(Mat mat,Vec v,PetscInt idx[])
4559: {
4566: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
4567: if (!mat->ops->getrowminabs) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
4568: MatCheckPreallocated(mat,1);
4569: if (idx) {PetscMemzero(idx,mat->rmap->n*sizeof(PetscInt));}
4571: (*mat->ops->getrowminabs)(mat,v,idx);
4572: PetscObjectStateIncrease((PetscObject)v);
4573: return(0);
4574: }
4576: /*@C
4577: MatGetRowMax - Gets the maximum value (of the real part) of each
4578: row of the matrix
4580: Logically Collective on Mat and Vec
4582: Input Parameters:
4583: . mat - the matrix
4585: Output Parameter:
4586: + v - the vector for storing the maximums
4587: - idx - the indices of the column found for each row (optional)
4589: Level: intermediate
4591: Notes: The result of this call are the same as if one converted the matrix to dense format
4592: and found the minimum value in each row (i.e. the implicit zeros are counted as zeros).
4594: This code is only implemented for a couple of matrix formats.
4596: Concepts: matrices^getting row maximums
4598: .seealso: MatGetDiagonal(), MatCreateSubMatrices(), MatCreateSubMatrix(), MatGetRowMaxAbs(), MatGetRowMin()
4599: @*/
4600: PetscErrorCode MatGetRowMax(Mat mat,Vec v,PetscInt idx[])
4601: {
4608: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
4609: if (!mat->ops->getrowmax) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
4610: MatCheckPreallocated(mat,1);
4612: (*mat->ops->getrowmax)(mat,v,idx);
4613: PetscObjectStateIncrease((PetscObject)v);
4614: return(0);
4615: }
4617: /*@C
4618: MatGetRowMaxAbs - Gets the maximum value (in absolute value) of each
4619: row of the matrix
4621: Logically Collective on Mat and Vec
4623: Input Parameters:
4624: . mat - the matrix
4626: Output Parameter:
4627: + v - the vector for storing the maximums
4628: - idx - the indices of the column found for each row (or NULL if not needed)
4630: Level: intermediate
4632: Notes: if a row is completely empty or has only 0.0 values then the idx[] value for that
4633: row is 0 (the first column).
4635: This code is only implemented for a couple of matrix formats.
4637: Concepts: matrices^getting row maximums
4639: .seealso: MatGetDiagonal(), MatCreateSubMatrices(), MatCreateSubMatrix(), MatGetRowMax(), MatGetRowMin()
4640: @*/
4641: PetscErrorCode MatGetRowMaxAbs(Mat mat,Vec v,PetscInt idx[])
4642: {
4649: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
4650: if (!mat->ops->getrowmaxabs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
4651: MatCheckPreallocated(mat,1);
4652: if (idx) {PetscMemzero(idx,mat->rmap->n*sizeof(PetscInt));}
4654: (*mat->ops->getrowmaxabs)(mat,v,idx);
4655: PetscObjectStateIncrease((PetscObject)v);
4656: return(0);
4657: }
4659: /*@
4660: MatGetRowSum - Gets the sum of each row of the matrix
4662: Logically or Neighborhood Collective on Mat and Vec
4664: Input Parameters:
4665: . mat - the matrix
4667: Output Parameter:
4668: . v - the vector for storing the sum of rows
4670: Level: intermediate
4672: Notes: This code is slow since it is not currently specialized for different formats
4674: Concepts: matrices^getting row sums
4676: .seealso: MatGetDiagonal(), MatCreateSubMatrices(), MatCreateSubMatrix(), MatGetRowMax(), MatGetRowMin()
4677: @*/
4678: PetscErrorCode MatGetRowSum(Mat mat, Vec v)
4679: {
4680: Vec ones;
4687: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
4688: MatCheckPreallocated(mat,1);
4689: MatCreateVecs(mat,&ones,NULL);
4690: VecSet(ones,1.);
4691: MatMult(mat,ones,v);
4692: VecDestroy(&ones);
4693: return(0);
4694: }
4696: /*@
4697: MatTranspose - Computes an in-place or out-of-place transpose of a matrix.
4699: Collective on Mat
4701: Input Parameter:
4702: + mat - the matrix to transpose
4703: - reuse - either MAT_INITIAL_MATRIX, MAT_REUSE_MATRIX, or MAT_INPLACE_MATRIX
4705: Output Parameters:
4706: . B - the transpose
4708: Notes:
4709: If you use MAT_INPLACE_MATRIX then you must pass in &mat for B
4711: MAT_REUSE_MATRIX causes the B matrix from a previous call to this function with MAT_INITIAL_MATRIX to be used
4713: Consider using MatCreateTranspose() instead if you only need a matrix that behaves like the transpose, but don't need the storage to be changed.
4715: Level: intermediate
4717: Concepts: matrices^transposing
4719: .seealso: MatMultTranspose(), MatMultTransposeAdd(), MatIsTranspose(), MatReuse
4720: @*/
4721: PetscErrorCode MatTranspose(Mat mat,MatReuse reuse,Mat *B)
4722: {
4728: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
4729: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
4730: if (!mat->ops->transpose) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
4731: if (reuse == MAT_INPLACE_MATRIX && mat != *B) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"MAT_INPLACE_MATRIX requires last matrix to match first");
4732: if (reuse == MAT_REUSE_MATRIX && mat == *B) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Perhaps you mean MAT_INPLACE_MATRIX");
4733: MatCheckPreallocated(mat,1);
4735: PetscLogEventBegin(MAT_Transpose,mat,0,0,0);
4736: (*mat->ops->transpose)(mat,reuse,B);
4737: PetscLogEventEnd(MAT_Transpose,mat,0,0,0);
4738: if (B) {PetscObjectStateIncrease((PetscObject)*B);}
4739: return(0);
4740: }
4742: /*@
4743: MatIsTranspose - Test whether a matrix is another one's transpose,
4744: or its own, in which case it tests symmetry.
4746: Collective on Mat
4748: Input Parameter:
4749: + A - the matrix to test
4750: - B - the matrix to test against, this can equal the first parameter
4752: Output Parameters:
4753: . flg - the result
4755: Notes:
4756: Only available for SeqAIJ/MPIAIJ matrices. The sequential algorithm
4757: has a running time of the order of the number of nonzeros; the parallel
4758: test involves parallel copies of the block-offdiagonal parts of the matrix.
4760: Level: intermediate
4762: Concepts: matrices^transposing, matrix^symmetry
4764: .seealso: MatTranspose(), MatIsSymmetric(), MatIsHermitian()
4765: @*/
4766: PetscErrorCode MatIsTranspose(Mat A,Mat B,PetscReal tol,PetscBool *flg)
4767: {
4768: PetscErrorCode ierr,(*f)(Mat,Mat,PetscReal,PetscBool*),(*g)(Mat,Mat,PetscReal,PetscBool*);
4774: PetscObjectQueryFunction((PetscObject)A,"MatIsTranspose_C",&f);
4775: PetscObjectQueryFunction((PetscObject)B,"MatIsTranspose_C",&g);
4776: *flg = PETSC_FALSE;
4777: if (f && g) {
4778: if (f == g) {
4779: (*f)(A,B,tol,flg);
4780: } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_NOTSAMETYPE,"Matrices do not have the same comparator for symmetry test");
4781: } else {
4782: MatType mattype;
4783: if (!f) {
4784: MatGetType(A,&mattype);
4785: } else {
4786: MatGetType(B,&mattype);
4787: }
4788: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix of type <%s> does not support checking for transpose",mattype);
4789: }
4790: return(0);
4791: }
4793: /*@
4794: MatHermitianTranspose - Computes an in-place or out-of-place transpose of a matrix in complex conjugate.
4796: Collective on Mat
4798: Input Parameter:
4799: + mat - the matrix to transpose and complex conjugate
4800: - reuse - MAT_INITIAL_MATRIX to create a new matrix, MAT_INPLACE_MATRIX to reuse the first argument to store the transpose
4802: Output Parameters:
4803: . B - the Hermitian
4805: Level: intermediate
4807: Concepts: matrices^transposing, complex conjugatex
4809: .seealso: MatTranspose(), MatMultTranspose(), MatMultTransposeAdd(), MatIsTranspose(), MatReuse
4810: @*/
4811: PetscErrorCode MatHermitianTranspose(Mat mat,MatReuse reuse,Mat *B)
4812: {
4816: MatTranspose(mat,reuse,B);
4817: #if defined(PETSC_USE_COMPLEX)
4818: MatConjugate(*B);
4819: #endif
4820: return(0);
4821: }
4823: /*@
4824: MatIsHermitianTranspose - Test whether a matrix is another one's Hermitian transpose,
4826: Collective on Mat
4828: Input Parameter:
4829: + A - the matrix to test
4830: - B - the matrix to test against, this can equal the first parameter
4832: Output Parameters:
4833: . flg - the result
4835: Notes:
4836: Only available for SeqAIJ/MPIAIJ matrices. The sequential algorithm
4837: has a running time of the order of the number of nonzeros; the parallel
4838: test involves parallel copies of the block-offdiagonal parts of the matrix.
4840: Level: intermediate
4842: Concepts: matrices^transposing, matrix^symmetry
4844: .seealso: MatTranspose(), MatIsSymmetric(), MatIsHermitian(), MatIsTranspose()
4845: @*/
4846: PetscErrorCode MatIsHermitianTranspose(Mat A,Mat B,PetscReal tol,PetscBool *flg)
4847: {
4848: PetscErrorCode ierr,(*f)(Mat,Mat,PetscReal,PetscBool*),(*g)(Mat,Mat,PetscReal,PetscBool*);
4854: PetscObjectQueryFunction((PetscObject)A,"MatIsHermitianTranspose_C",&f);
4855: PetscObjectQueryFunction((PetscObject)B,"MatIsHermitianTranspose_C",&g);
4856: if (f && g) {
4857: if (f==g) {
4858: (*f)(A,B,tol,flg);
4859: } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_NOTSAMETYPE,"Matrices do not have the same comparator for Hermitian test");
4860: }
4861: return(0);
4862: }
4864: /*@
4865: MatPermute - Creates a new matrix with rows and columns permuted from the
4866: original.
4868: Collective on Mat
4870: Input Parameters:
4871: + mat - the matrix to permute
4872: . row - row permutation, each processor supplies only the permutation for its rows
4873: - col - column permutation, each processor supplies only the permutation for its columns
4875: Output Parameters:
4876: . B - the permuted matrix
4878: Level: advanced
4880: Note:
4881: The index sets map from row/col of permuted matrix to row/col of original matrix.
4882: The index sets should be on the same communicator as Mat and have the same local sizes.
4884: Concepts: matrices^permuting
4886: .seealso: MatGetOrdering(), ISAllGather()
4888: @*/
4889: PetscErrorCode MatPermute(Mat mat,IS row,IS col,Mat *B)
4890: {
4899: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
4900: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
4901: if (!mat->ops->permute) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatPermute not available for Mat type %s",((PetscObject)mat)->type_name);
4902: MatCheckPreallocated(mat,1);
4904: (*mat->ops->permute)(mat,row,col,B);
4905: PetscObjectStateIncrease((PetscObject)*B);
4906: return(0);
4907: }
4909: /*@
4910: MatEqual - Compares two matrices.
4912: Collective on Mat
4914: Input Parameters:
4915: + A - the first matrix
4916: - B - the second matrix
4918: Output Parameter:
4919: . flg - PETSC_TRUE if the matrices are equal; PETSC_FALSE otherwise.
4921: Level: intermediate
4923: Concepts: matrices^equality between
4924: @*/
4925: PetscErrorCode MatEqual(Mat A,Mat B,PetscBool *flg)
4926: {
4936: MatCheckPreallocated(B,2);
4937: if (!A->assembled) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
4938: if (!B->assembled) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
4939: if (A->rmap->N != B->rmap->N || A->cmap->N != B->cmap->N) SETERRQ4(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Mat A,Mat B: global dim %D %D %D %D",A->rmap->N,B->rmap->N,A->cmap->N,B->cmap->N);
4940: if (!A->ops->equal) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Mat type %s",((PetscObject)A)->type_name);
4941: if (!B->ops->equal) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Mat type %s",((PetscObject)B)->type_name);
4942: if (A->ops->equal != B->ops->equal) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"A is type: %s\nB is type: %s",((PetscObject)A)->type_name,((PetscObject)B)->type_name);
4943: MatCheckPreallocated(A,1);
4945: (*A->ops->equal)(A,B,flg);
4946: return(0);
4947: }
4949: /*@C
4950: MatDiagonalScale - Scales a matrix on the left and right by diagonal
4951: matrices that are stored as vectors. Either of the two scaling
4952: matrices can be NULL.
4954: Collective on Mat
4956: Input Parameters:
4957: + mat - the matrix to be scaled
4958: . l - the left scaling vector (or NULL)
4959: - r - the right scaling vector (or NULL)
4961: Notes:
4962: MatDiagonalScale() computes A = LAR, where
4963: L = a diagonal matrix (stored as a vector), R = a diagonal matrix (stored as a vector)
4964: The L scales the rows of the matrix, the R scales the columns of the matrix.
4966: Level: intermediate
4968: Concepts: matrices^diagonal scaling
4969: Concepts: diagonal scaling of matrices
4971: .seealso: MatScale(), MatShift(), MatDiagonalSet()
4972: @*/
4973: PetscErrorCode MatDiagonalScale(Mat mat,Vec l,Vec r)
4974: {
4980: if (!mat->ops->diagonalscale) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
4983: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
4984: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
4985: MatCheckPreallocated(mat,1);
4987: PetscLogEventBegin(MAT_Scale,mat,0,0,0);
4988: (*mat->ops->diagonalscale)(mat,l,r);
4989: PetscLogEventEnd(MAT_Scale,mat,0,0,0);
4990: PetscObjectStateIncrease((PetscObject)mat);
4991: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_VECCUDA)
4992: if (mat->valid_GPU_matrix != PETSC_OFFLOAD_UNALLOCATED) {
4993: mat->valid_GPU_matrix = PETSC_OFFLOAD_CPU;
4994: }
4995: #endif
4996: return(0);
4997: }
4999: /*@
5000: MatScale - Scales all elements of a matrix by a given number.
5002: Logically Collective on Mat
5004: Input Parameters:
5005: + mat - the matrix to be scaled
5006: - a - the scaling value
5008: Output Parameter:
5009: . mat - the scaled matrix
5011: Level: intermediate
5013: Concepts: matrices^scaling all entries
5015: .seealso: MatDiagonalScale()
5016: @*/
5017: PetscErrorCode MatScale(Mat mat,PetscScalar a)
5018: {
5024: if (a != (PetscScalar)1.0 && !mat->ops->scale) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
5025: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
5026: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
5028: MatCheckPreallocated(mat,1);
5030: PetscLogEventBegin(MAT_Scale,mat,0,0,0);
5031: if (a != (PetscScalar)1.0) {
5032: (*mat->ops->scale)(mat,a);
5033: PetscObjectStateIncrease((PetscObject)mat);
5034: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_VECCUDA)
5035: if (mat->valid_GPU_matrix != PETSC_OFFLOAD_UNALLOCATED) {
5036: mat->valid_GPU_matrix = PETSC_OFFLOAD_CPU;
5037: }
5038: #endif
5039: }
5040: PetscLogEventEnd(MAT_Scale,mat,0,0,0);
5041: return(0);
5042: }
5044: static PetscErrorCode MatNorm_Basic(Mat A,NormType type,PetscReal *nrm)
5045: {
5049: if (type == NORM_1 || type == NORM_INFINITY) {
5050: Vec l,r;
5052: MatCreateVecs(A,&r,&l);
5053: if (type == NORM_INFINITY) {
5054: VecSet(r,1.);
5055: MatMult(A,r,l);
5056: VecNorm(l,NORM_INFINITY,nrm);
5057: } else {
5058: VecSet(l,1.);
5059: MatMultTranspose(A,l,r);
5060: VecNorm(r,NORM_INFINITY,nrm);
5061: }
5062: VecDestroy(&l);
5063: VecDestroy(&r);
5064: } else SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Matrix class %s, norm type %d",((PetscObject)A)->type_name,type);
5065: return(0);
5066: }
5068: /*@
5069: MatNorm - Calculates various norms of a matrix.
5071: Collective on Mat
5073: Input Parameters:
5074: + mat - the matrix
5075: - type - the type of norm, NORM_1, NORM_FROBENIUS, NORM_INFINITY
5077: Output Parameters:
5078: . nrm - the resulting norm
5080: Level: intermediate
5082: Concepts: matrices^norm
5083: Concepts: norm^of matrix
5084: @*/
5085: PetscErrorCode MatNorm(Mat mat,NormType type,PetscReal *nrm)
5086: {
5095: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
5096: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
5097: MatCheckPreallocated(mat,1);
5099: if (!mat->ops->norm) {
5100: MatNorm_Basic(mat,type,nrm);
5101: } else {
5102: (*mat->ops->norm)(mat,type,nrm);
5103: }
5104: return(0);
5105: }
5107: /*
5108: This variable is used to prevent counting of MatAssemblyBegin() that
5109: are called from within a MatAssemblyEnd().
5110: */
5111: static PetscInt MatAssemblyEnd_InUse = 0;
5112: /*@
5113: MatAssemblyBegin - Begins assembling the matrix. This routine should
5114: be called after completing all calls to MatSetValues().
5116: Collective on Mat
5118: Input Parameters:
5119: + mat - the matrix
5120: - type - type of assembly, either MAT_FLUSH_ASSEMBLY or MAT_FINAL_ASSEMBLY
5122: Notes:
5123: MatSetValues() generally caches the values. The matrix is ready to
5124: use only after MatAssemblyBegin() and MatAssemblyEnd() have been called.
5125: Use MAT_FLUSH_ASSEMBLY when switching between ADD_VALUES and INSERT_VALUES
5126: in MatSetValues(); use MAT_FINAL_ASSEMBLY for the final assembly before
5127: using the matrix.
5129: ALL processes that share a matrix MUST call MatAssemblyBegin() and MatAssemblyEnd() the SAME NUMBER of times, and each time with the
5130: same flag of MAT_FLUSH_ASSEMBLY or MAT_FINAL_ASSEMBLY for all processes. Thus you CANNOT locally change from ADD_VALUES to INSERT_VALUES, that is
5131: a global collective operation requring all processes that share the matrix.
5133: Space for preallocated nonzeros that is not filled by a call to MatSetValues() or a related routine are compressed
5134: out by assembly. If you intend to use that extra space on a subsequent assembly, be sure to insert explicit zeros
5135: before MAT_FINAL_ASSEMBLY so the space is not compressed out.
5137: Level: beginner
5139: Concepts: matrices^assembling
5141: .seealso: MatAssemblyEnd(), MatSetValues(), MatAssembled()
5142: @*/
5143: PetscErrorCode MatAssemblyBegin(Mat mat,MatAssemblyType type)
5144: {
5150: MatCheckPreallocated(mat,1);
5151: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix.\nDid you forget to call MatSetUnfactored()?");
5152: if (mat->assembled) {
5153: mat->was_assembled = PETSC_TRUE;
5154: mat->assembled = PETSC_FALSE;
5155: }
5156: if (!MatAssemblyEnd_InUse) {
5157: PetscLogEventBegin(MAT_AssemblyBegin,mat,0,0,0);
5158: if (mat->ops->assemblybegin) {(*mat->ops->assemblybegin)(mat,type);}
5159: PetscLogEventEnd(MAT_AssemblyBegin,mat,0,0,0);
5160: } else if (mat->ops->assemblybegin) {
5161: (*mat->ops->assemblybegin)(mat,type);
5162: }
5163: return(0);
5164: }
5166: /*@
5167: MatAssembled - Indicates if a matrix has been assembled and is ready for
5168: use; for example, in matrix-vector product.
5170: Not Collective
5172: Input Parameter:
5173: . mat - the matrix
5175: Output Parameter:
5176: . assembled - PETSC_TRUE or PETSC_FALSE
5178: Level: advanced
5180: Concepts: matrices^assembled?
5182: .seealso: MatAssemblyEnd(), MatSetValues(), MatAssemblyBegin()
5183: @*/
5184: PetscErrorCode MatAssembled(Mat mat,PetscBool *assembled)
5185: {
5190: *assembled = mat->assembled;
5191: return(0);
5192: }
5194: /*@
5195: MatAssemblyEnd - Completes assembling the matrix. This routine should
5196: be called after MatAssemblyBegin().
5198: Collective on Mat
5200: Input Parameters:
5201: + mat - the matrix
5202: - type - type of assembly, either MAT_FLUSH_ASSEMBLY or MAT_FINAL_ASSEMBLY
5204: Options Database Keys:
5205: + -mat_view ::ascii_info - Prints info on matrix at conclusion of MatEndAssembly()
5206: . -mat_view ::ascii_info_detail - Prints more detailed info
5207: . -mat_view - Prints matrix in ASCII format
5208: . -mat_view ::ascii_matlab - Prints matrix in Matlab format
5209: . -mat_view draw - PetscDraws nonzero structure of matrix, using MatView() and PetscDrawOpenX().
5210: . -display <name> - Sets display name (default is host)
5211: . -draw_pause <sec> - Sets number of seconds to pause after display
5212: . -mat_view socket - Sends matrix to socket, can be accessed from Matlab (See Users-Manual: Chapter 12 Using MATLAB with PETSc )
5213: . -viewer_socket_machine <machine> - Machine to use for socket
5214: . -viewer_socket_port <port> - Port number to use for socket
5215: - -mat_view binary:filename[:append] - Save matrix to file in binary format
5217: Notes:
5218: MatSetValues() generally caches the values. The matrix is ready to
5219: use only after MatAssemblyBegin() and MatAssemblyEnd() have been called.
5220: Use MAT_FLUSH_ASSEMBLY when switching between ADD_VALUES and INSERT_VALUES
5221: in MatSetValues(); use MAT_FINAL_ASSEMBLY for the final assembly before
5222: using the matrix.
5224: Space for preallocated nonzeros that is not filled by a call to MatSetValues() or a related routine are compressed
5225: out by assembly. If you intend to use that extra space on a subsequent assembly, be sure to insert explicit zeros
5226: before MAT_FINAL_ASSEMBLY so the space is not compressed out.
5228: Level: beginner
5230: .seealso: MatAssemblyBegin(), MatSetValues(), PetscDrawOpenX(), PetscDrawCreate(), MatView(), MatAssembled(), PetscViewerSocketOpen()
5231: @*/
5232: PetscErrorCode MatAssemblyEnd(Mat mat,MatAssemblyType type)
5233: {
5234: PetscErrorCode ierr;
5235: static PetscInt inassm = 0;
5236: PetscBool flg = PETSC_FALSE;
5242: inassm++;
5243: MatAssemblyEnd_InUse++;
5244: if (MatAssemblyEnd_InUse == 1) { /* Do the logging only the first time through */
5245: PetscLogEventBegin(MAT_AssemblyEnd,mat,0,0,0);
5246: if (mat->ops->assemblyend) {
5247: (*mat->ops->assemblyend)(mat,type);
5248: }
5249: PetscLogEventEnd(MAT_AssemblyEnd,mat,0,0,0);
5250: } else if (mat->ops->assemblyend) {
5251: (*mat->ops->assemblyend)(mat,type);
5252: }
5254: /* Flush assembly is not a true assembly */
5255: if (type != MAT_FLUSH_ASSEMBLY) {
5256: mat->assembled = PETSC_TRUE; mat->num_ass++;
5257: }
5258: mat->insertmode = NOT_SET_VALUES;
5259: MatAssemblyEnd_InUse--;
5260: PetscObjectStateIncrease((PetscObject)mat);
5261: if (!mat->symmetric_eternal) {
5262: mat->symmetric_set = PETSC_FALSE;
5263: mat->hermitian_set = PETSC_FALSE;
5264: mat->structurally_symmetric_set = PETSC_FALSE;
5265: }
5266: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_VECCUDA)
5267: if (mat->valid_GPU_matrix != PETSC_OFFLOAD_UNALLOCATED) {
5268: mat->valid_GPU_matrix = PETSC_OFFLOAD_CPU;
5269: }
5270: #endif
5271: if (inassm == 1 && type != MAT_FLUSH_ASSEMBLY) {
5272: MatViewFromOptions(mat,NULL,"-mat_view");
5274: if (mat->checksymmetryonassembly) {
5275: MatIsSymmetric(mat,mat->checksymmetrytol,&flg);
5276: if (flg) {
5277: PetscPrintf(PetscObjectComm((PetscObject)mat),"Matrix is symmetric (tolerance %g)\n",(double)mat->checksymmetrytol);
5278: } else {
5279: PetscPrintf(PetscObjectComm((PetscObject)mat),"Matrix is not symmetric (tolerance %g)\n",(double)mat->checksymmetrytol);
5280: }
5281: }
5282: if (mat->nullsp && mat->checknullspaceonassembly) {
5283: MatNullSpaceTest(mat->nullsp,mat,NULL);
5284: }
5285: }
5286: inassm--;
5287: return(0);
5288: }
5290: /*@
5291: MatSetOption - Sets a parameter option for a matrix. Some options
5292: may be specific to certain storage formats. Some options
5293: determine how values will be inserted (or added). Sorted,
5294: row-oriented input will generally assemble the fastest. The default
5295: is row-oriented.
5297: Logically Collective on Mat for certain operations, such as MAT_SPD, not collective for MAT_ROW_ORIENTED, see MatOption
5299: Input Parameters:
5300: + mat - the matrix
5301: . option - the option, one of those listed below (and possibly others),
5302: - flg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE)
5304: Options Describing Matrix Structure:
5305: + MAT_SPD - symmetric positive definite
5306: . MAT_SYMMETRIC - symmetric in terms of both structure and value
5307: . MAT_HERMITIAN - transpose is the complex conjugation
5308: . MAT_STRUCTURALLY_SYMMETRIC - symmetric nonzero structure
5309: - MAT_SYMMETRY_ETERNAL - if you would like the symmetry/Hermitian flag
5310: you set to be kept with all future use of the matrix
5311: including after MatAssemblyBegin/End() which could
5312: potentially change the symmetry structure, i.e. you
5313: KNOW the matrix will ALWAYS have the property you set.
5316: Options For Use with MatSetValues():
5317: Insert a logically dense subblock, which can be
5318: . MAT_ROW_ORIENTED - row-oriented (default)
5320: Note these options reflect the data you pass in with MatSetValues(); it has
5321: nothing to do with how the data is stored internally in the matrix
5322: data structure.
5324: When (re)assembling a matrix, we can restrict the input for
5325: efficiency/debugging purposes. These options include:
5326: + MAT_NEW_NONZERO_LOCATIONS - additional insertions will be allowed if they generate a new nonzero (slow)
5327: . MAT_NEW_DIAGONALS - new diagonals will be allowed (for block diagonal format only)
5328: . MAT_IGNORE_OFF_PROC_ENTRIES - drops off-processor entries
5329: . MAT_NEW_NONZERO_LOCATION_ERR - generates an error for new matrix entry
5330: . MAT_USE_HASH_TABLE - uses a hash table to speed up matrix assembly
5331: . MAT_NO_OFF_PROC_ENTRIES - you know each process will only set values for its own rows, will generate an error if
5332: any process sets values for another process. This avoids all reductions in the MatAssembly routines and thus improves
5333: performance for very large process counts.
5334: - MAT_SUBSET_OFF_PROC_ENTRIES - you know that the first assembly after setting this flag will set a superset
5335: of the off-process entries required for all subsequent assemblies. This avoids a rendezvous step in the MatAssembly
5336: functions, instead sending only neighbor messages.
5338: Notes:
5339: Except for MAT_UNUSED_NONZERO_LOCATION_ERR and MAT_ROW_ORIENTED all processes that share the matrix must pass the same value in flg!
5341: Some options are relevant only for particular matrix types and
5342: are thus ignored by others. Other options are not supported by
5343: certain matrix types and will generate an error message if set.
5345: If using a Fortran 77 module to compute a matrix, one may need to
5346: use the column-oriented option (or convert to the row-oriented
5347: format).
5349: MAT_NEW_NONZERO_LOCATIONS set to PETSC_FALSE indicates that any add or insertion
5350: that would generate a new entry in the nonzero structure is instead
5351: ignored. Thus, if memory has not alredy been allocated for this particular
5352: data, then the insertion is ignored. For dense matrices, in which
5353: the entire array is allocated, no entries are ever ignored.
5354: Set after the first MatAssemblyEnd(). If this option is set then the MatAssemblyBegin/End() processes has one less global reduction
5356: MAT_NEW_NONZERO_LOCATION_ERR set to PETSC_TRUE indicates that any add or insertion
5357: that would generate a new entry in the nonzero structure instead produces
5358: an error. (Currently supported for AIJ and BAIJ formats only.) If this option is set then the MatAssemblyBegin/End() processes has one less global reduction
5360: MAT_NEW_NONZERO_ALLOCATION_ERR set to PETSC_TRUE indicates that any add or insertion
5361: that would generate a new entry that has not been preallocated will
5362: instead produce an error. (Currently supported for AIJ and BAIJ formats
5363: only.) This is a useful flag when debugging matrix memory preallocation.
5364: If this option is set then the MatAssemblyBegin/End() processes has one less global reduction
5366: MAT_IGNORE_OFF_PROC_ENTRIES set to PETSC_TRUE indicates entries destined for
5367: other processors should be dropped, rather than stashed.
5368: This is useful if you know that the "owning" processor is also
5369: always generating the correct matrix entries, so that PETSc need
5370: not transfer duplicate entries generated on another processor.
5372: MAT_USE_HASH_TABLE indicates that a hash table be used to improve the
5373: searches during matrix assembly. When this flag is set, the hash table
5374: is created during the first Matrix Assembly. This hash table is
5375: used the next time through, during MatSetVaules()/MatSetVaulesBlocked()
5376: to improve the searching of indices. MAT_NEW_NONZERO_LOCATIONS flag
5377: should be used with MAT_USE_HASH_TABLE flag. This option is currently
5378: supported by MATMPIBAIJ format only.
5380: MAT_KEEP_NONZERO_PATTERN indicates when MatZeroRows() is called the zeroed entries
5381: are kept in the nonzero structure
5383: MAT_IGNORE_ZERO_ENTRIES - for AIJ/IS matrices this will stop zero values from creating
5384: a zero location in the matrix
5386: MAT_USE_INODES - indicates using inode version of the code - works with AIJ matrix types
5388: MAT_NO_OFF_PROC_ZERO_ROWS - you know each process will only zero its own rows. This avoids all reductions in the
5389: zero row routines and thus improves performance for very large process counts.
5391: MAT_IGNORE_LOWER_TRIANGULAR - For SBAIJ matrices will ignore any insertions you make in the lower triangular
5392: part of the matrix (since they should match the upper triangular part).
5394: Notes: Can only be called after MatSetSizes() and MatSetType() have been set.
5396: Level: intermediate
5398: Concepts: matrices^setting options
5400: .seealso: MatOption, Mat
5402: @*/
5403: PetscErrorCode MatSetOption(Mat mat,MatOption op,PetscBool flg)
5404: {
5410: if (op > 0) {
5413: }
5415: if (((int) op) <= MAT_OPTION_MIN || ((int) op) >= MAT_OPTION_MAX) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_OUTOFRANGE,"Options %d is out of range",(int)op);
5416: if (!((PetscObject)mat)->type_name) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_TYPENOTSET,"Cannot set options until type and size have been set, see MatSetType() and MatSetSizes()");
5418: switch (op) {
5419: case MAT_NO_OFF_PROC_ENTRIES:
5420: mat->nooffprocentries = flg;
5421: return(0);
5422: break;
5423: case MAT_SUBSET_OFF_PROC_ENTRIES:
5424: mat->subsetoffprocentries = flg;
5425: return(0);
5426: case MAT_NO_OFF_PROC_ZERO_ROWS:
5427: mat->nooffproczerorows = flg;
5428: return(0);
5429: break;
5430: case MAT_SPD:
5431: mat->spd_set = PETSC_TRUE;
5432: mat->spd = flg;
5433: if (flg) {
5434: mat->symmetric = PETSC_TRUE;
5435: mat->structurally_symmetric = PETSC_TRUE;
5436: mat->symmetric_set = PETSC_TRUE;
5437: mat->structurally_symmetric_set = PETSC_TRUE;
5438: }
5439: break;
5440: case MAT_SYMMETRIC:
5441: mat->symmetric = flg;
5442: if (flg) mat->structurally_symmetric = PETSC_TRUE;
5443: mat->symmetric_set = PETSC_TRUE;
5444: mat->structurally_symmetric_set = flg;
5445: #if !defined(PETSC_USE_COMPLEX)
5446: mat->hermitian = flg;
5447: mat->hermitian_set = PETSC_TRUE;
5448: #endif
5449: break;
5450: case MAT_HERMITIAN:
5451: mat->hermitian = flg;
5452: if (flg) mat->structurally_symmetric = PETSC_TRUE;
5453: mat->hermitian_set = PETSC_TRUE;
5454: mat->structurally_symmetric_set = flg;
5455: #if !defined(PETSC_USE_COMPLEX)
5456: mat->symmetric = flg;
5457: mat->symmetric_set = PETSC_TRUE;
5458: #endif
5459: break;
5460: case MAT_STRUCTURALLY_SYMMETRIC:
5461: mat->structurally_symmetric = flg;
5462: mat->structurally_symmetric_set = PETSC_TRUE;
5463: break;
5464: case MAT_SYMMETRY_ETERNAL:
5465: mat->symmetric_eternal = flg;
5466: break;
5467: case MAT_STRUCTURE_ONLY:
5468: mat->structure_only = flg;
5469: break;
5470: default:
5471: break;
5472: }
5473: if (mat->ops->setoption) {
5474: (*mat->ops->setoption)(mat,op,flg);
5475: }
5476: return(0);
5477: }
5479: /*@
5480: MatGetOption - Gets a parameter option that has been set for a matrix.
5482: Logically Collective on Mat for certain operations, such as MAT_SPD, not collective for MAT_ROW_ORIENTED, see MatOption
5484: Input Parameters:
5485: + mat - the matrix
5486: - option - the option, this only responds to certain options, check the code for which ones
5488: Output Parameter:
5489: . flg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE)
5491: Notes: Can only be called after MatSetSizes() and MatSetType() have been set.
5493: Level: intermediate
5495: Concepts: matrices^setting options
5497: .seealso: MatOption, MatSetOption()
5499: @*/
5500: PetscErrorCode MatGetOption(Mat mat,MatOption op,PetscBool *flg)
5501: {
5506: if (((int) op) <= MAT_OPTION_MIN || ((int) op) >= MAT_OPTION_MAX) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_OUTOFRANGE,"Options %d is out of range",(int)op);
5507: if (!((PetscObject)mat)->type_name) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_TYPENOTSET,"Cannot get options until type and size have been set, see MatSetType() and MatSetSizes()");
5509: switch (op) {
5510: case MAT_NO_OFF_PROC_ENTRIES:
5511: *flg = mat->nooffprocentries;
5512: break;
5513: case MAT_NO_OFF_PROC_ZERO_ROWS:
5514: *flg = mat->nooffproczerorows;
5515: break;
5516: case MAT_SYMMETRIC:
5517: *flg = mat->symmetric;
5518: break;
5519: case MAT_HERMITIAN:
5520: *flg = mat->hermitian;
5521: break;
5522: case MAT_STRUCTURALLY_SYMMETRIC:
5523: *flg = mat->structurally_symmetric;
5524: break;
5525: case MAT_SYMMETRY_ETERNAL:
5526: *flg = mat->symmetric_eternal;
5527: break;
5528: case MAT_SPD:
5529: *flg = mat->spd;
5530: break;
5531: default:
5532: break;
5533: }
5534: return(0);
5535: }
5537: /*@
5538: MatZeroEntries - Zeros all entries of a matrix. For sparse matrices
5539: this routine retains the old nonzero structure.
5541: Logically Collective on Mat
5543: Input Parameters:
5544: . mat - the matrix
5546: Level: intermediate
5548: Notes: If the matrix was not preallocated then a default, likely poor preallocation will be set in the matrix, so this should be called after the preallocation phase.
5549: See the Performance chapter of the users manual for information on preallocating matrices.
5551: Concepts: matrices^zeroing
5553: .seealso: MatZeroRows()
5554: @*/
5555: PetscErrorCode MatZeroEntries(Mat mat)
5556: {
5562: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
5563: if (mat->insertmode != NOT_SET_VALUES) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for matrices where you have set values but not yet assembled");
5564: if (!mat->ops->zeroentries) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
5565: MatCheckPreallocated(mat,1);
5567: PetscLogEventBegin(MAT_ZeroEntries,mat,0,0,0);
5568: (*mat->ops->zeroentries)(mat);
5569: PetscLogEventEnd(MAT_ZeroEntries,mat,0,0,0);
5570: PetscObjectStateIncrease((PetscObject)mat);
5571: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_VECCUDA)
5572: if (mat->valid_GPU_matrix != PETSC_OFFLOAD_UNALLOCATED) {
5573: mat->valid_GPU_matrix = PETSC_OFFLOAD_CPU;
5574: }
5575: #endif
5576: return(0);
5577: }
5579: /*@C
5580: MatZeroRowsColumns - Zeros all entries (except possibly the main diagonal)
5581: of a set of rows and columns of a matrix.
5583: Collective on Mat
5585: Input Parameters:
5586: + mat - the matrix
5587: . numRows - the number of rows to remove
5588: . rows - the global row indices
5589: . diag - value put in all diagonals of eliminated rows (0.0 will even eliminate diagonal entry)
5590: . x - optional vector of solutions for zeroed rows (other entries in vector are not used)
5591: - b - optional vector of right hand side, that will be adjusted by provided solution
5593: Notes:
5594: This does not change the nonzero structure of the matrix, it merely zeros those entries in the matrix.
5596: The user can set a value in the diagonal entry (or for the AIJ and
5597: row formats can optionally remove the main diagonal entry from the
5598: nonzero structure as well, by passing 0.0 as the final argument).
5600: For the parallel case, all processes that share the matrix (i.e.,
5601: those in the communicator used for matrix creation) MUST call this
5602: routine, regardless of whether any rows being zeroed are owned by
5603: them.
5605: Each processor can indicate any rows in the entire matrix to be zeroed (i.e. each process does NOT have to
5606: list only rows local to itself).
5608: The option MAT_NO_OFF_PROC_ZERO_ROWS does not apply to this routine.
5610: Level: intermediate
5612: Concepts: matrices^zeroing rows
5614: .seealso: MatZeroRowsIS(), MatZeroRows(), MatZeroRowsLocalIS(), MatZeroRowsStencil(), MatZeroEntries(), MatZeroRowsLocal(), MatSetOption(),
5615: MatZeroRowsColumnsLocal(), MatZeroRowsColumnsLocalIS(), MatZeroRowsColumnsIS(), MatZeroRowsColumnsStencil()
5616: @*/
5617: PetscErrorCode MatZeroRowsColumns(Mat mat,PetscInt numRows,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
5618: {
5625: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
5626: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
5627: if (!mat->ops->zerorowscolumns) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
5628: MatCheckPreallocated(mat,1);
5630: (*mat->ops->zerorowscolumns)(mat,numRows,rows,diag,x,b);
5631: MatViewFromOptions(mat,NULL,"-mat_view");
5632: PetscObjectStateIncrease((PetscObject)mat);
5633: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_VECCUDA)
5634: if (mat->valid_GPU_matrix != PETSC_OFFLOAD_UNALLOCATED) {
5635: mat->valid_GPU_matrix = PETSC_OFFLOAD_CPU;
5636: }
5637: #endif
5638: return(0);
5639: }
5641: /*@C
5642: MatZeroRowsColumnsIS - Zeros all entries (except possibly the main diagonal)
5643: of a set of rows and columns of a matrix.
5645: Collective on Mat
5647: Input Parameters:
5648: + mat - the matrix
5649: . is - the rows to zero
5650: . diag - value put in all diagonals of eliminated rows (0.0 will even eliminate diagonal entry)
5651: . x - optional vector of solutions for zeroed rows (other entries in vector are not used)
5652: - b - optional vector of right hand side, that will be adjusted by provided solution
5654: Notes:
5655: This does not change the nonzero structure of the matrix, it merely zeros those entries in the matrix.
5657: The user can set a value in the diagonal entry (or for the AIJ and
5658: row formats can optionally remove the main diagonal entry from the
5659: nonzero structure as well, by passing 0.0 as the final argument).
5661: For the parallel case, all processes that share the matrix (i.e.,
5662: those in the communicator used for matrix creation) MUST call this
5663: routine, regardless of whether any rows being zeroed are owned by
5664: them.
5666: Each processor can indicate any rows in the entire matrix to be zeroed (i.e. each process does NOT have to
5667: list only rows local to itself).
5669: The option MAT_NO_OFF_PROC_ZERO_ROWS does not apply to this routine.
5671: Level: intermediate
5673: Concepts: matrices^zeroing rows
5675: .seealso: MatZeroRowsIS(), MatZeroRowsColumns(), MatZeroRowsLocalIS(), MatZeroRowsStencil(), MatZeroEntries(), MatZeroRowsLocal(), MatSetOption(),
5676: MatZeroRowsColumnsLocal(), MatZeroRowsColumnsLocalIS(), MatZeroRows(), MatZeroRowsColumnsStencil()
5677: @*/
5678: PetscErrorCode MatZeroRowsColumnsIS(Mat mat,IS is,PetscScalar diag,Vec x,Vec b)
5679: {
5681: PetscInt numRows;
5682: const PetscInt *rows;
5689: ISGetLocalSize(is,&numRows);
5690: ISGetIndices(is,&rows);
5691: MatZeroRowsColumns(mat,numRows,rows,diag,x,b);
5692: ISRestoreIndices(is,&rows);
5693: return(0);
5694: }
5696: /*@C
5697: MatZeroRows - Zeros all entries (except possibly the main diagonal)
5698: of a set of rows of a matrix.
5700: Collective on Mat
5702: Input Parameters:
5703: + mat - the matrix
5704: . numRows - the number of rows to remove
5705: . rows - the global row indices
5706: . diag - value put in all diagonals of eliminated rows (0.0 will even eliminate diagonal entry)
5707: . x - optional vector of solutions for zeroed rows (other entries in vector are not used)
5708: - b - optional vector of right hand side, that will be adjusted by provided solution
5710: Notes:
5711: For the AIJ and BAIJ matrix formats this removes the old nonzero structure,
5712: but does not release memory. For the dense and block diagonal
5713: formats this does not alter the nonzero structure.
5715: If the option MatSetOption(mat,MAT_KEEP_NONZERO_PATTERN,PETSC_TRUE) the nonzero structure
5716: of the matrix is not changed (even for AIJ and BAIJ matrices) the values are
5717: merely zeroed.
5719: The user can set a value in the diagonal entry (or for the AIJ and
5720: row formats can optionally remove the main diagonal entry from the
5721: nonzero structure as well, by passing 0.0 as the final argument).
5723: For the parallel case, all processes that share the matrix (i.e.,
5724: those in the communicator used for matrix creation) MUST call this
5725: routine, regardless of whether any rows being zeroed are owned by
5726: them.
5728: Each processor can indicate any rows in the entire matrix to be zeroed (i.e. each process does NOT have to
5729: list only rows local to itself).
5731: You can call MatSetOption(mat,MAT_NO_OFF_PROC_ZERO_ROWS,PETSC_TRUE) if each process indicates only rows it
5732: owns that are to be zeroed. This saves a global synchronization in the implementation.
5734: Level: intermediate
5736: Concepts: matrices^zeroing rows
5738: .seealso: MatZeroRowsIS(), MatZeroRowsColumns(), MatZeroRowsLocalIS(), MatZeroRowsStencil(), MatZeroEntries(), MatZeroRowsLocal(), MatSetOption(),
5739: MatZeroRowsColumnsLocal(), MatZeroRowsColumnsLocalIS(), MatZeroRowsColumnsIS(), MatZeroRowsColumnsStencil()
5740: @*/
5741: PetscErrorCode MatZeroRows(Mat mat,PetscInt numRows,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
5742: {
5749: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
5750: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
5751: if (!mat->ops->zerorows) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
5752: MatCheckPreallocated(mat,1);
5754: (*mat->ops->zerorows)(mat,numRows,rows,diag,x,b);
5755: MatViewFromOptions(mat,NULL,"-mat_view");
5756: PetscObjectStateIncrease((PetscObject)mat);
5757: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_VECCUDA)
5758: if (mat->valid_GPU_matrix != PETSC_OFFLOAD_UNALLOCATED) {
5759: mat->valid_GPU_matrix = PETSC_OFFLOAD_CPU;
5760: }
5761: #endif
5762: return(0);
5763: }
5765: /*@C
5766: MatZeroRowsIS - Zeros all entries (except possibly the main diagonal)
5767: of a set of rows of a matrix.
5769: Collective on Mat
5771: Input Parameters:
5772: + mat - the matrix
5773: . is - index set of rows to remove
5774: . diag - value put in all diagonals of eliminated rows
5775: . x - optional vector of solutions for zeroed rows (other entries in vector are not used)
5776: - b - optional vector of right hand side, that will be adjusted by provided solution
5778: Notes:
5779: For the AIJ and BAIJ matrix formats this removes the old nonzero structure,
5780: but does not release memory. For the dense and block diagonal
5781: formats this does not alter the nonzero structure.
5783: If the option MatSetOption(mat,MAT_KEEP_NONZERO_PATTERN,PETSC_TRUE) the nonzero structure
5784: of the matrix is not changed (even for AIJ and BAIJ matrices) the values are
5785: merely zeroed.
5787: The user can set a value in the diagonal entry (or for the AIJ and
5788: row formats can optionally remove the main diagonal entry from the
5789: nonzero structure as well, by passing 0.0 as the final argument).
5791: For the parallel case, all processes that share the matrix (i.e.,
5792: those in the communicator used for matrix creation) MUST call this
5793: routine, regardless of whether any rows being zeroed are owned by
5794: them.
5796: Each processor can indicate any rows in the entire matrix to be zeroed (i.e. each process does NOT have to
5797: list only rows local to itself).
5799: You can call MatSetOption(mat,MAT_NO_OFF_PROC_ZERO_ROWS,PETSC_TRUE) if each process indicates only rows it
5800: owns that are to be zeroed. This saves a global synchronization in the implementation.
5802: Level: intermediate
5804: Concepts: matrices^zeroing rows
5806: .seealso: MatZeroRows(), MatZeroRowsColumns(), MatZeroRowsLocalIS(), MatZeroRowsStencil(), MatZeroEntries(), MatZeroRowsLocal(), MatSetOption(),
5807: MatZeroRowsColumnsLocal(), MatZeroRowsColumnsLocalIS(), MatZeroRowsColumnsIS(), MatZeroRowsColumnsStencil()
5808: @*/
5809: PetscErrorCode MatZeroRowsIS(Mat mat,IS is,PetscScalar diag,Vec x,Vec b)
5810: {
5811: PetscInt numRows;
5812: const PetscInt *rows;
5819: ISGetLocalSize(is,&numRows);
5820: ISGetIndices(is,&rows);
5821: MatZeroRows(mat,numRows,rows,diag,x,b);
5822: ISRestoreIndices(is,&rows);
5823: return(0);
5824: }
5826: /*@C
5827: MatZeroRowsStencil - Zeros all entries (except possibly the main diagonal)
5828: of a set of rows of a matrix. These rows must be local to the process.
5830: Collective on Mat
5832: Input Parameters:
5833: + mat - the matrix
5834: . numRows - the number of rows to remove
5835: . rows - the grid coordinates (and component number when dof > 1) for matrix rows
5836: . diag - value put in all diagonals of eliminated rows (0.0 will even eliminate diagonal entry)
5837: . x - optional vector of solutions for zeroed rows (other entries in vector are not used)
5838: - b - optional vector of right hand side, that will be adjusted by provided solution
5840: Notes:
5841: For the AIJ and BAIJ matrix formats this removes the old nonzero structure,
5842: but does not release memory. For the dense and block diagonal
5843: formats this does not alter the nonzero structure.
5845: If the option MatSetOption(mat,MAT_KEEP_NONZERO_PATTERN,PETSC_TRUE) the nonzero structure
5846: of the matrix is not changed (even for AIJ and BAIJ matrices) the values are
5847: merely zeroed.
5849: The user can set a value in the diagonal entry (or for the AIJ and
5850: row formats can optionally remove the main diagonal entry from the
5851: nonzero structure as well, by passing 0.0 as the final argument).
5853: For the parallel case, all processes that share the matrix (i.e.,
5854: those in the communicator used for matrix creation) MUST call this
5855: routine, regardless of whether any rows being zeroed are owned by
5856: them.
5858: Each processor can indicate any rows in the entire matrix to be zeroed (i.e. each process does NOT have to
5859: list only rows local to itself).
5861: The grid coordinates are across the entire grid, not just the local portion
5863: In Fortran idxm and idxn should be declared as
5864: $ MatStencil idxm(4,m)
5865: and the values inserted using
5866: $ idxm(MatStencil_i,1) = i
5867: $ idxm(MatStencil_j,1) = j
5868: $ idxm(MatStencil_k,1) = k
5869: $ idxm(MatStencil_c,1) = c
5870: etc
5872: For periodic boundary conditions use negative indices for values to the left (below 0; that are to be
5873: obtained by wrapping values from right edge). For values to the right of the last entry using that index plus one
5874: etc to obtain values that obtained by wrapping the values from the left edge. This does not work for anything but the
5875: DM_BOUNDARY_PERIODIC boundary type.
5877: For indices that don't mean anything for your case (like the k index when working in 2d) or the c index when you have
5878: a single value per point) you can skip filling those indices.
5880: Level: intermediate
5882: Concepts: matrices^zeroing rows
5884: .seealso: MatZeroRowsIS(), MatZeroRowsColumns(), MatZeroRowsLocalIS(), MatZeroRowsl(), MatZeroEntries(), MatZeroRowsLocal(), MatSetOption(),
5885: MatZeroRowsColumnsLocal(), MatZeroRowsColumnsLocalIS(), MatZeroRowsColumnsIS(), MatZeroRowsColumnsStencil()
5886: @*/
5887: PetscErrorCode MatZeroRowsStencil(Mat mat,PetscInt numRows,const MatStencil rows[],PetscScalar diag,Vec x,Vec b)
5888: {
5889: PetscInt dim = mat->stencil.dim;
5890: PetscInt sdim = dim - (1 - (PetscInt) mat->stencil.noc);
5891: PetscInt *dims = mat->stencil.dims+1;
5892: PetscInt *starts = mat->stencil.starts;
5893: PetscInt *dxm = (PetscInt*) rows;
5894: PetscInt *jdxm, i, j, tmp, numNewRows = 0;
5902: PetscMalloc1(numRows, &jdxm);
5903: for (i = 0; i < numRows; ++i) {
5904: /* Skip unused dimensions (they are ordered k, j, i, c) */
5905: for (j = 0; j < 3-sdim; ++j) dxm++;
5906: /* Local index in X dir */
5907: tmp = *dxm++ - starts[0];
5908: /* Loop over remaining dimensions */
5909: for (j = 0; j < dim-1; ++j) {
5910: /* If nonlocal, set index to be negative */
5911: if ((*dxm++ - starts[j+1]) < 0 || tmp < 0) tmp = PETSC_MIN_INT;
5912: /* Update local index */
5913: else tmp = tmp*dims[j] + *(dxm-1) - starts[j+1];
5914: }
5915: /* Skip component slot if necessary */
5916: if (mat->stencil.noc) dxm++;
5917: /* Local row number */
5918: if (tmp >= 0) {
5919: jdxm[numNewRows++] = tmp;
5920: }
5921: }
5922: MatZeroRowsLocal(mat,numNewRows,jdxm,diag,x,b);
5923: PetscFree(jdxm);
5924: return(0);
5925: }
5927: /*@C
5928: MatZeroRowsColumnsStencil - Zeros all row and column entries (except possibly the main diagonal)
5929: of a set of rows and columns of a matrix.
5931: Collective on Mat
5933: Input Parameters:
5934: + mat - the matrix
5935: . numRows - the number of rows/columns to remove
5936: . rows - the grid coordinates (and component number when dof > 1) for matrix rows
5937: . diag - value put in all diagonals of eliminated rows (0.0 will even eliminate diagonal entry)
5938: . x - optional vector of solutions for zeroed rows (other entries in vector are not used)
5939: - b - optional vector of right hand side, that will be adjusted by provided solution
5941: Notes:
5942: For the AIJ and BAIJ matrix formats this removes the old nonzero structure,
5943: but does not release memory. For the dense and block diagonal
5944: formats this does not alter the nonzero structure.
5946: If the option MatSetOption(mat,MAT_KEEP_NONZERO_PATTERN,PETSC_TRUE) the nonzero structure
5947: of the matrix is not changed (even for AIJ and BAIJ matrices) the values are
5948: merely zeroed.
5950: The user can set a value in the diagonal entry (or for the AIJ and
5951: row formats can optionally remove the main diagonal entry from the
5952: nonzero structure as well, by passing 0.0 as the final argument).
5954: For the parallel case, all processes that share the matrix (i.e.,
5955: those in the communicator used for matrix creation) MUST call this
5956: routine, regardless of whether any rows being zeroed are owned by
5957: them.
5959: Each processor can indicate any rows in the entire matrix to be zeroed (i.e. each process does NOT have to
5960: list only rows local to itself, but the row/column numbers are given in local numbering).
5962: The grid coordinates are across the entire grid, not just the local portion
5964: In Fortran idxm and idxn should be declared as
5965: $ MatStencil idxm(4,m)
5966: and the values inserted using
5967: $ idxm(MatStencil_i,1) = i
5968: $ idxm(MatStencil_j,1) = j
5969: $ idxm(MatStencil_k,1) = k
5970: $ idxm(MatStencil_c,1) = c
5971: etc
5973: For periodic boundary conditions use negative indices for values to the left (below 0; that are to be
5974: obtained by wrapping values from right edge). For values to the right of the last entry using that index plus one
5975: etc to obtain values that obtained by wrapping the values from the left edge. This does not work for anything but the
5976: DM_BOUNDARY_PERIODIC boundary type.
5978: For indices that don't mean anything for your case (like the k index when working in 2d) or the c index when you have
5979: a single value per point) you can skip filling those indices.
5981: Level: intermediate
5983: Concepts: matrices^zeroing rows
5985: .seealso: MatZeroRowsIS(), MatZeroRowsColumns(), MatZeroRowsLocalIS(), MatZeroRowsStencil(), MatZeroEntries(), MatZeroRowsLocal(), MatSetOption(),
5986: MatZeroRowsColumnsLocal(), MatZeroRowsColumnsLocalIS(), MatZeroRowsColumnsIS(), MatZeroRows()
5987: @*/
5988: PetscErrorCode MatZeroRowsColumnsStencil(Mat mat,PetscInt numRows,const MatStencil rows[],PetscScalar diag,Vec x,Vec b)
5989: {
5990: PetscInt dim = mat->stencil.dim;
5991: PetscInt sdim = dim - (1 - (PetscInt) mat->stencil.noc);
5992: PetscInt *dims = mat->stencil.dims+1;
5993: PetscInt *starts = mat->stencil.starts;
5994: PetscInt *dxm = (PetscInt*) rows;
5995: PetscInt *jdxm, i, j, tmp, numNewRows = 0;
6003: PetscMalloc1(numRows, &jdxm);
6004: for (i = 0; i < numRows; ++i) {
6005: /* Skip unused dimensions (they are ordered k, j, i, c) */
6006: for (j = 0; j < 3-sdim; ++j) dxm++;
6007: /* Local index in X dir */
6008: tmp = *dxm++ - starts[0];
6009: /* Loop over remaining dimensions */
6010: for (j = 0; j < dim-1; ++j) {
6011: /* If nonlocal, set index to be negative */
6012: if ((*dxm++ - starts[j+1]) < 0 || tmp < 0) tmp = PETSC_MIN_INT;
6013: /* Update local index */
6014: else tmp = tmp*dims[j] + *(dxm-1) - starts[j+1];
6015: }
6016: /* Skip component slot if necessary */
6017: if (mat->stencil.noc) dxm++;
6018: /* Local row number */
6019: if (tmp >= 0) {
6020: jdxm[numNewRows++] = tmp;
6021: }
6022: }
6023: MatZeroRowsColumnsLocal(mat,numNewRows,jdxm,diag,x,b);
6024: PetscFree(jdxm);
6025: return(0);
6026: }
6028: /*@C
6029: MatZeroRowsLocal - Zeros all entries (except possibly the main diagonal)
6030: of a set of rows of a matrix; using local numbering of rows.
6032: Collective on Mat
6034: Input Parameters:
6035: + mat - the matrix
6036: . numRows - the number of rows to remove
6037: . rows - the global row indices
6038: . diag - value put in all diagonals of eliminated rows
6039: . x - optional vector of solutions for zeroed rows (other entries in vector are not used)
6040: - b - optional vector of right hand side, that will be adjusted by provided solution
6042: Notes:
6043: Before calling MatZeroRowsLocal(), the user must first set the
6044: local-to-global mapping by calling MatSetLocalToGlobalMapping().
6046: For the AIJ matrix formats this removes the old nonzero structure,
6047: but does not release memory. For the dense and block diagonal
6048: formats this does not alter the nonzero structure.
6050: If the option MatSetOption(mat,MAT_KEEP_NONZERO_PATTERN,PETSC_TRUE) the nonzero structure
6051: of the matrix is not changed (even for AIJ and BAIJ matrices) the values are
6052: merely zeroed.
6054: The user can set a value in the diagonal entry (or for the AIJ and
6055: row formats can optionally remove the main diagonal entry from the
6056: nonzero structure as well, by passing 0.0 as the final argument).
6058: You can call MatSetOption(mat,MAT_NO_OFF_PROC_ZERO_ROWS,PETSC_TRUE) if each process indicates only rows it
6059: owns that are to be zeroed. This saves a global synchronization in the implementation.
6061: Level: intermediate
6063: Concepts: matrices^zeroing
6065: .seealso: MatZeroRowsIS(), MatZeroRowsColumns(), MatZeroRowsLocalIS(), MatZeroRowsStencil(), MatZeroEntries(), MatZeroRows(), MatSetOption(),
6066: MatZeroRowsColumnsLocal(), MatZeroRowsColumnsLocalIS(), MatZeroRowsColumnsIS(), MatZeroRowsColumnsStencil()
6067: @*/
6068: PetscErrorCode MatZeroRowsLocal(Mat mat,PetscInt numRows,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
6069: {
6076: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
6077: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
6078: MatCheckPreallocated(mat,1);
6080: if (mat->ops->zerorowslocal) {
6081: (*mat->ops->zerorowslocal)(mat,numRows,rows,diag,x,b);
6082: } else {
6083: IS is, newis;
6084: const PetscInt *newRows;
6086: if (!mat->rmap->mapping) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Need to provide local to global mapping to matrix first");
6087: ISCreateGeneral(PETSC_COMM_SELF,numRows,rows,PETSC_COPY_VALUES,&is);
6088: ISLocalToGlobalMappingApplyIS(mat->rmap->mapping,is,&newis);
6089: ISGetIndices(newis,&newRows);
6090: (*mat->ops->zerorows)(mat,numRows,newRows,diag,x,b);
6091: ISRestoreIndices(newis,&newRows);
6092: ISDestroy(&newis);
6093: ISDestroy(&is);
6094: }
6095: PetscObjectStateIncrease((PetscObject)mat);
6096: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_VECCUDA)
6097: if (mat->valid_GPU_matrix != PETSC_OFFLOAD_UNALLOCATED) {
6098: mat->valid_GPU_matrix = PETSC_OFFLOAD_CPU;
6099: }
6100: #endif
6101: return(0);
6102: }
6104: /*@C
6105: MatZeroRowsLocalIS - Zeros all entries (except possibly the main diagonal)
6106: of a set of rows of a matrix; using local numbering of rows.
6108: Collective on Mat
6110: Input Parameters:
6111: + mat - the matrix
6112: . is - index set of rows to remove
6113: . diag - value put in all diagonals of eliminated rows
6114: . x - optional vector of solutions for zeroed rows (other entries in vector are not used)
6115: - b - optional vector of right hand side, that will be adjusted by provided solution
6117: Notes:
6118: Before calling MatZeroRowsLocalIS(), the user must first set the
6119: local-to-global mapping by calling MatSetLocalToGlobalMapping().
6121: For the AIJ matrix formats this removes the old nonzero structure,
6122: but does not release memory. For the dense and block diagonal
6123: formats this does not alter the nonzero structure.
6125: If the option MatSetOption(mat,MAT_KEEP_NONZERO_PATTERN,PETSC_TRUE) the nonzero structure
6126: of the matrix is not changed (even for AIJ and BAIJ matrices) the values are
6127: merely zeroed.
6129: The user can set a value in the diagonal entry (or for the AIJ and
6130: row formats can optionally remove the main diagonal entry from the
6131: nonzero structure as well, by passing 0.0 as the final argument).
6133: You can call MatSetOption(mat,MAT_NO_OFF_PROC_ZERO_ROWS,PETSC_TRUE) if each process indicates only rows it
6134: owns that are to be zeroed. This saves a global synchronization in the implementation.
6136: Level: intermediate
6138: Concepts: matrices^zeroing
6140: .seealso: MatZeroRowsIS(), MatZeroRowsColumns(), MatZeroRows(), MatZeroRowsStencil(), MatZeroEntries(), MatZeroRowsLocal(), MatSetOption(),
6141: MatZeroRowsColumnsLocal(), MatZeroRowsColumnsLocalIS(), MatZeroRowsColumnsIS(), MatZeroRowsColumnsStencil()
6142: @*/
6143: PetscErrorCode MatZeroRowsLocalIS(Mat mat,IS is,PetscScalar diag,Vec x,Vec b)
6144: {
6146: PetscInt numRows;
6147: const PetscInt *rows;
6153: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
6154: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
6155: MatCheckPreallocated(mat,1);
6157: ISGetLocalSize(is,&numRows);
6158: ISGetIndices(is,&rows);
6159: MatZeroRowsLocal(mat,numRows,rows,diag,x,b);
6160: ISRestoreIndices(is,&rows);
6161: return(0);
6162: }
6164: /*@C
6165: MatZeroRowsColumnsLocal - Zeros all entries (except possibly the main diagonal)
6166: of a set of rows and columns of a matrix; using local numbering of rows.
6168: Collective on Mat
6170: Input Parameters:
6171: + mat - the matrix
6172: . numRows - the number of rows to remove
6173: . rows - the global row indices
6174: . diag - value put in all diagonals of eliminated rows
6175: . x - optional vector of solutions for zeroed rows (other entries in vector are not used)
6176: - b - optional vector of right hand side, that will be adjusted by provided solution
6178: Notes:
6179: Before calling MatZeroRowsColumnsLocal(), the user must first set the
6180: local-to-global mapping by calling MatSetLocalToGlobalMapping().
6182: The user can set a value in the diagonal entry (or for the AIJ and
6183: row formats can optionally remove the main diagonal entry from the
6184: nonzero structure as well, by passing 0.0 as the final argument).
6186: Level: intermediate
6188: Concepts: matrices^zeroing
6190: .seealso: MatZeroRowsIS(), MatZeroRowsColumns(), MatZeroRowsLocalIS(), MatZeroRowsStencil(), MatZeroEntries(), MatZeroRowsLocal(), MatSetOption(),
6191: MatZeroRows(), MatZeroRowsColumnsLocalIS(), MatZeroRowsColumnsIS(), MatZeroRowsColumnsStencil()
6192: @*/
6193: PetscErrorCode MatZeroRowsColumnsLocal(Mat mat,PetscInt numRows,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
6194: {
6196: IS is, newis;
6197: const PetscInt *newRows;
6203: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
6204: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
6205: MatCheckPreallocated(mat,1);
6207: if (!mat->cmap->mapping) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Need to provide local to global mapping to matrix first");
6208: ISCreateGeneral(PETSC_COMM_SELF,numRows,rows,PETSC_COPY_VALUES,&is);
6209: ISLocalToGlobalMappingApplyIS(mat->cmap->mapping,is,&newis);
6210: ISGetIndices(newis,&newRows);
6211: (*mat->ops->zerorowscolumns)(mat,numRows,newRows,diag,x,b);
6212: ISRestoreIndices(newis,&newRows);
6213: ISDestroy(&newis);
6214: ISDestroy(&is);
6215: PetscObjectStateIncrease((PetscObject)mat);
6216: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_VECCUDA)
6217: if (mat->valid_GPU_matrix != PETSC_OFFLOAD_UNALLOCATED) {
6218: mat->valid_GPU_matrix = PETSC_OFFLOAD_CPU;
6219: }
6220: #endif
6221: return(0);
6222: }
6224: /*@C
6225: MatZeroRowsColumnsLocalIS - Zeros all entries (except possibly the main diagonal)
6226: of a set of rows and columns of a matrix; using local numbering of rows.
6228: Collective on Mat
6230: Input Parameters:
6231: + mat - the matrix
6232: . is - index set of rows to remove
6233: . diag - value put in all diagonals of eliminated rows
6234: . x - optional vector of solutions for zeroed rows (other entries in vector are not used)
6235: - b - optional vector of right hand side, that will be adjusted by provided solution
6237: Notes:
6238: Before calling MatZeroRowsColumnsLocalIS(), the user must first set the
6239: local-to-global mapping by calling MatSetLocalToGlobalMapping().
6241: The user can set a value in the diagonal entry (or for the AIJ and
6242: row formats can optionally remove the main diagonal entry from the
6243: nonzero structure as well, by passing 0.0 as the final argument).
6245: Level: intermediate
6247: Concepts: matrices^zeroing
6249: .seealso: MatZeroRowsIS(), MatZeroRowsColumns(), MatZeroRowsLocalIS(), MatZeroRowsStencil(), MatZeroEntries(), MatZeroRowsLocal(), MatSetOption(),
6250: MatZeroRowsColumnsLocal(), MatZeroRows(), MatZeroRowsColumnsIS(), MatZeroRowsColumnsStencil()
6251: @*/
6252: PetscErrorCode MatZeroRowsColumnsLocalIS(Mat mat,IS is,PetscScalar diag,Vec x,Vec b)
6253: {
6255: PetscInt numRows;
6256: const PetscInt *rows;
6262: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
6263: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
6264: MatCheckPreallocated(mat,1);
6266: ISGetLocalSize(is,&numRows);
6267: ISGetIndices(is,&rows);
6268: MatZeroRowsColumnsLocal(mat,numRows,rows,diag,x,b);
6269: ISRestoreIndices(is,&rows);
6270: return(0);
6271: }
6273: /*@C
6274: MatGetSize - Returns the numbers of rows and columns in a matrix.
6276: Not Collective
6278: Input Parameter:
6279: . mat - the matrix
6281: Output Parameters:
6282: + m - the number of global rows
6283: - n - the number of global columns
6285: Note: both output parameters can be NULL on input.
6287: Level: beginner
6289: Concepts: matrices^size
6291: .seealso: MatGetLocalSize()
6292: @*/
6293: PetscErrorCode MatGetSize(Mat mat,PetscInt *m,PetscInt *n)
6294: {
6297: if (m) *m = mat->rmap->N;
6298: if (n) *n = mat->cmap->N;
6299: return(0);
6300: }
6302: /*@C
6303: MatGetLocalSize - Returns the number of rows and columns in a matrix
6304: stored locally. This information may be implementation dependent, so
6305: use with care.
6307: Not Collective
6309: Input Parameters:
6310: . mat - the matrix
6312: Output Parameters:
6313: + m - the number of local rows
6314: - n - the number of local columns
6316: Note: both output parameters can be NULL on input.
6318: Level: beginner
6320: Concepts: matrices^local size
6322: .seealso: MatGetSize()
6323: @*/
6324: PetscErrorCode MatGetLocalSize(Mat mat,PetscInt *m,PetscInt *n)
6325: {
6330: if (m) *m = mat->rmap->n;
6331: if (n) *n = mat->cmap->n;
6332: return(0);
6333: }
6335: /*@C
6336: MatGetOwnershipRangeColumn - Returns the range of matrix columns associated with rows of a vector one multiplies by that owned by
6337: this processor. (The columns of the "diagonal block")
6339: Not Collective, unless matrix has not been allocated, then collective on Mat
6341: Input Parameters:
6342: . mat - the matrix
6344: Output Parameters:
6345: + m - the global index of the first local column
6346: - n - one more than the global index of the last local column
6348: Notes: both output parameters can be NULL on input.
6350: Level: developer
6352: Concepts: matrices^column ownership
6354: .seealso: MatGetOwnershipRange(), MatGetOwnershipRanges(), MatGetOwnershipRangesColumn()
6356: @*/
6357: PetscErrorCode MatGetOwnershipRangeColumn(Mat mat,PetscInt *m,PetscInt *n)
6358: {
6364: MatCheckPreallocated(mat,1);
6365: if (m) *m = mat->cmap->rstart;
6366: if (n) *n = mat->cmap->rend;
6367: return(0);
6368: }
6370: /*@C
6371: MatGetOwnershipRange - Returns the range of matrix rows owned by
6372: this processor, assuming that the matrix is laid out with the first
6373: n1 rows on the first processor, the next n2 rows on the second, etc.
6374: For certain parallel layouts this range may not be well defined.
6376: Not Collective
6378: Input Parameters:
6379: . mat - the matrix
6381: Output Parameters:
6382: + m - the global index of the first local row
6383: - n - one more than the global index of the last local row
6385: Note: Both output parameters can be NULL on input.
6386: $ This function requires that the matrix be preallocated. If you have not preallocated, consider using
6387: $ PetscSplitOwnership(MPI_Comm comm, PetscInt *n, PetscInt *N)
6388: $ and then MPI_Scan() to calculate prefix sums of the local sizes.
6390: Level: beginner
6392: Concepts: matrices^row ownership
6394: .seealso: MatGetOwnershipRanges(), MatGetOwnershipRangeColumn(), MatGetOwnershipRangesColumn(), PetscSplitOwnership(), PetscSplitOwnershipBlock()
6396: @*/
6397: PetscErrorCode MatGetOwnershipRange(Mat mat,PetscInt *m,PetscInt *n)
6398: {
6404: MatCheckPreallocated(mat,1);
6405: if (m) *m = mat->rmap->rstart;
6406: if (n) *n = mat->rmap->rend;
6407: return(0);
6408: }
6410: /*@C
6411: MatGetOwnershipRanges - Returns the range of matrix rows owned by
6412: each process
6414: Not Collective, unless matrix has not been allocated, then collective on Mat
6416: Input Parameters:
6417: . mat - the matrix
6419: Output Parameters:
6420: . ranges - start of each processors portion plus one more than the total length at the end
6422: Level: beginner
6424: Concepts: matrices^row ownership
6426: .seealso: MatGetOwnershipRange(), MatGetOwnershipRangeColumn(), MatGetOwnershipRangesColumn()
6428: @*/
6429: PetscErrorCode MatGetOwnershipRanges(Mat mat,const PetscInt **ranges)
6430: {
6436: MatCheckPreallocated(mat,1);
6437: PetscLayoutGetRanges(mat->rmap,ranges);
6438: return(0);
6439: }
6441: /*@C
6442: MatGetOwnershipRangesColumn - Returns the range of matrix columns associated with rows of a vector one multiplies by that owned by
6443: this processor. (The columns of the "diagonal blocks" for each process)
6445: Not Collective, unless matrix has not been allocated, then collective on Mat
6447: Input Parameters:
6448: . mat - the matrix
6450: Output Parameters:
6451: . ranges - start of each processors portion plus one more then the total length at the end
6453: Level: beginner
6455: Concepts: matrices^column ownership
6457: .seealso: MatGetOwnershipRange(), MatGetOwnershipRangeColumn(), MatGetOwnershipRanges()
6459: @*/
6460: PetscErrorCode MatGetOwnershipRangesColumn(Mat mat,const PetscInt **ranges)
6461: {
6467: MatCheckPreallocated(mat,1);
6468: PetscLayoutGetRanges(mat->cmap,ranges);
6469: return(0);
6470: }
6472: /*@C
6473: MatGetOwnershipIS - Get row and column ownership as index sets
6475: Not Collective
6477: Input Arguments:
6478: . A - matrix of type Elemental
6480: Output Arguments:
6481: + rows - rows in which this process owns elements
6482: . cols - columns in which this process owns elements
6484: Level: intermediate
6486: .seealso: MatGetOwnershipRange(), MatGetOwnershipRangeColumn(), MatSetValues(), MATELEMENTAL
6487: @*/
6488: PetscErrorCode MatGetOwnershipIS(Mat A,IS *rows,IS *cols)
6489: {
6490: PetscErrorCode ierr,(*f)(Mat,IS*,IS*);
6493: MatCheckPreallocated(A,1);
6494: PetscObjectQueryFunction((PetscObject)A,"MatGetOwnershipIS_C",&f);
6495: if (f) {
6496: (*f)(A,rows,cols);
6497: } else { /* Create a standard row-based partition, each process is responsible for ALL columns in their row block */
6498: if (rows) {ISCreateStride(PETSC_COMM_SELF,A->rmap->n,A->rmap->rstart,1,rows);}
6499: if (cols) {ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,cols);}
6500: }
6501: return(0);
6502: }
6504: /*@C
6505: MatILUFactorSymbolic - Performs symbolic ILU factorization of a matrix.
6506: Uses levels of fill only, not drop tolerance. Use MatLUFactorNumeric()
6507: to complete the factorization.
6509: Collective on Mat
6511: Input Parameters:
6512: + mat - the matrix
6513: . row - row permutation
6514: . column - column permutation
6515: - info - structure containing
6516: $ levels - number of levels of fill.
6517: $ expected fill - as ratio of original fill.
6518: $ 1 or 0 - indicating force fill on diagonal (improves robustness for matrices
6519: missing diagonal entries)
6521: Output Parameters:
6522: . fact - new matrix that has been symbolically factored
6524: Notes: See Users-Manual: ch_mat for additional information about choosing the fill factor for better efficiency.
6526: Most users should employ the simplified KSP interface for linear solvers
6527: instead of working directly with matrix algebra routines such as this.
6528: See, e.g., KSPCreate().
6530: Level: developer
6532: Concepts: matrices^symbolic LU factorization
6533: Concepts: matrices^factorization
6534: Concepts: LU^symbolic factorization
6536: .seealso: MatLUFactorSymbolic(), MatLUFactorNumeric(), MatCholeskyFactor()
6537: MatGetOrdering(), MatFactorInfo
6539: Developer Note: fortran interface is not autogenerated as the f90
6540: interface defintion cannot be generated correctly [due to MatFactorInfo]
6542: @*/
6543: PetscErrorCode MatILUFactorSymbolic(Mat fact,Mat mat,IS row,IS col,const MatFactorInfo *info)
6544: {
6554: if (info->levels < 0) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_OUTOFRANGE,"Levels of fill negative %D",(PetscInt)info->levels);
6555: if (info->fill < 1.0) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_OUTOFRANGE,"Expected fill less than 1.0 %g",(double)info->fill);
6556: if (!(fact)->ops->ilufactorsymbolic) {
6557: MatSolverType spackage;
6558: MatFactorGetSolverType(fact,&spackage);
6559: SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Matrix type %s symbolic ILU using solver package %s",((PetscObject)mat)->type_name,spackage);
6560: }
6561: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
6562: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
6563: MatCheckPreallocated(mat,2);
6565: PetscLogEventBegin(MAT_ILUFactorSymbolic,mat,row,col,0);
6566: (fact->ops->ilufactorsymbolic)(fact,mat,row,col,info);
6567: PetscLogEventEnd(MAT_ILUFactorSymbolic,mat,row,col,0);
6568: return(0);
6569: }
6571: /*@C
6572: MatICCFactorSymbolic - Performs symbolic incomplete
6573: Cholesky factorization for a symmetric matrix. Use
6574: MatCholeskyFactorNumeric() to complete the factorization.
6576: Collective on Mat
6578: Input Parameters:
6579: + mat - the matrix
6580: . perm - row and column permutation
6581: - info - structure containing
6582: $ levels - number of levels of fill.
6583: $ expected fill - as ratio of original fill.
6585: Output Parameter:
6586: . fact - the factored matrix
6588: Notes:
6589: Most users should employ the KSP interface for linear solvers
6590: instead of working directly with matrix algebra routines such as this.
6591: See, e.g., KSPCreate().
6593: Level: developer
6595: Concepts: matrices^symbolic incomplete Cholesky factorization
6596: Concepts: matrices^factorization
6597: Concepts: Cholsky^symbolic factorization
6599: .seealso: MatCholeskyFactorNumeric(), MatCholeskyFactor(), MatFactorInfo
6601: Developer Note: fortran interface is not autogenerated as the f90
6602: interface defintion cannot be generated correctly [due to MatFactorInfo]
6604: @*/
6605: PetscErrorCode MatICCFactorSymbolic(Mat fact,Mat mat,IS perm,const MatFactorInfo *info)
6606: {
6615: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
6616: if (info->levels < 0) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_OUTOFRANGE,"Levels negative %D",(PetscInt) info->levels);
6617: if (info->fill < 1.0) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_OUTOFRANGE,"Expected fill less than 1.0 %g",(double)info->fill);
6618: if (!(fact)->ops->iccfactorsymbolic) {
6619: MatSolverType spackage;
6620: MatFactorGetSolverType(fact,&spackage);
6621: SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Matrix type %s symbolic ICC using solver package %s",((PetscObject)mat)->type_name,spackage);
6622: }
6623: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
6624: MatCheckPreallocated(mat,2);
6626: PetscLogEventBegin(MAT_ICCFactorSymbolic,mat,perm,0,0);
6627: (fact->ops->iccfactorsymbolic)(fact,mat,perm,info);
6628: PetscLogEventEnd(MAT_ICCFactorSymbolic,mat,perm,0,0);
6629: return(0);
6630: }
6632: /*@C
6633: MatCreateSubMatrices - Extracts several submatrices from a matrix. If submat
6634: points to an array of valid matrices, they may be reused to store the new
6635: submatrices.
6637: Collective on Mat
6639: Input Parameters:
6640: + mat - the matrix
6641: . n - the number of submatrixes to be extracted (on this processor, may be zero)
6642: . irow, icol - index sets of rows and columns to extract
6643: - scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
6645: Output Parameter:
6646: . submat - the array of submatrices
6648: Notes:
6649: MatCreateSubMatrices() can extract ONLY sequential submatrices
6650: (from both sequential and parallel matrices). Use MatCreateSubMatrix()
6651: to extract a parallel submatrix.
6653: Some matrix types place restrictions on the row and column
6654: indices, such as that they be sorted or that they be equal to each other.
6656: The index sets may not have duplicate entries.
6658: When extracting submatrices from a parallel matrix, each processor can
6659: form a different submatrix by setting the rows and columns of its
6660: individual index sets according to the local submatrix desired.
6662: When finished using the submatrices, the user should destroy
6663: them with MatDestroySubMatrices().
6665: MAT_REUSE_MATRIX can only be used when the nonzero structure of the
6666: original matrix has not changed from that last call to MatCreateSubMatrices().
6668: This routine creates the matrices in submat; you should NOT create them before
6669: calling it. It also allocates the array of matrix pointers submat.
6671: For BAIJ matrices the index sets must respect the block structure, that is if they
6672: request one row/column in a block, they must request all rows/columns that are in
6673: that block. For example, if the block size is 2 you cannot request just row 0 and
6674: column 0.
6676: Fortran Note:
6677: The Fortran interface is slightly different from that given below; it
6678: requires one to pass in as submat a Mat (integer) array of size at least n+1.
6680: Level: advanced
6682: Concepts: matrices^accessing submatrices
6683: Concepts: submatrices
6685: .seealso: MatDestroySubMatrices(), MatCreateSubMatrix(), MatGetRow(), MatGetDiagonal(), MatReuse
6686: @*/
6687: PetscErrorCode MatCreateSubMatrices(Mat mat,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *submat[])
6688: {
6690: PetscInt i;
6691: PetscBool eq;
6696: if (n) {
6701: }
6703: if (n && scall == MAT_REUSE_MATRIX) {
6706: }
6707: if (!mat->ops->createsubmatrices) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
6708: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
6709: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
6710: MatCheckPreallocated(mat,1);
6712: PetscLogEventBegin(MAT_CreateSubMats,mat,0,0,0);
6713: (*mat->ops->createsubmatrices)(mat,n,irow,icol,scall,submat);
6714: PetscLogEventEnd(MAT_CreateSubMats,mat,0,0,0);
6715: for (i=0; i<n; i++) {
6716: (*submat)[i]->factortype = MAT_FACTOR_NONE; /* in case in place factorization was previously done on submatrix */
6717: if (mat->symmetric || mat->structurally_symmetric || mat->hermitian) {
6718: ISEqual(irow[i],icol[i],&eq);
6719: if (eq) {
6720: if (mat->symmetric) {
6721: MatSetOption((*submat)[i],MAT_SYMMETRIC,PETSC_TRUE);
6722: } else if (mat->hermitian) {
6723: MatSetOption((*submat)[i],MAT_HERMITIAN,PETSC_TRUE);
6724: } else if (mat->structurally_symmetric) {
6725: MatSetOption((*submat)[i],MAT_STRUCTURALLY_SYMMETRIC,PETSC_TRUE);
6726: }
6727: }
6728: }
6729: }
6730: return(0);
6731: }
6733: /*@C
6734: MatCreateSubMatricesMPI - Extracts MPI submatrices across a sub communicator of mat (by pairs of IS that may live on subcomms).
6736: Collective on Mat
6738: Input Parameters:
6739: + mat - the matrix
6740: . n - the number of submatrixes to be extracted
6741: . irow, icol - index sets of rows and columns to extract
6742: - scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
6744: Output Parameter:
6745: . submat - the array of submatrices
6747: Level: advanced
6749: Concepts: matrices^accessing submatrices
6750: Concepts: submatrices
6752: .seealso: MatCreateSubMatrices(), MatCreateSubMatrix(), MatGetRow(), MatGetDiagonal(), MatReuse
6753: @*/
6754: PetscErrorCode MatCreateSubMatricesMPI(Mat mat,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *submat[])
6755: {
6757: PetscInt i;
6758: PetscBool eq;
6763: if (n) {
6768: }
6770: if (n && scall == MAT_REUSE_MATRIX) {
6773: }
6774: if (!mat->ops->createsubmatricesmpi) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
6775: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
6776: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
6777: MatCheckPreallocated(mat,1);
6779: PetscLogEventBegin(MAT_CreateSubMats,mat,0,0,0);
6780: (*mat->ops->createsubmatricesmpi)(mat,n,irow,icol,scall,submat);
6781: PetscLogEventEnd(MAT_CreateSubMats,mat,0,0,0);
6782: for (i=0; i<n; i++) {
6783: if (mat->symmetric || mat->structurally_symmetric || mat->hermitian) {
6784: ISEqual(irow[i],icol[i],&eq);
6785: if (eq) {
6786: if (mat->symmetric) {
6787: MatSetOption((*submat)[i],MAT_SYMMETRIC,PETSC_TRUE);
6788: } else if (mat->hermitian) {
6789: MatSetOption((*submat)[i],MAT_HERMITIAN,PETSC_TRUE);
6790: } else if (mat->structurally_symmetric) {
6791: MatSetOption((*submat)[i],MAT_STRUCTURALLY_SYMMETRIC,PETSC_TRUE);
6792: }
6793: }
6794: }
6795: }
6796: return(0);
6797: }
6799: /*@C
6800: MatDestroyMatrices - Destroys an array of matrices.
6802: Collective on Mat
6804: Input Parameters:
6805: + n - the number of local matrices
6806: - mat - the matrices (note that this is a pointer to the array of matrices)
6808: Level: advanced
6810: Notes: Frees not only the matrices, but also the array that contains the matrices
6811: In Fortran will not free the array.
6813: .seealso: MatCreateSubMatrices() MatDestroySubMatrices()
6814: @*/
6815: PetscErrorCode MatDestroyMatrices(PetscInt n,Mat *mat[])
6816: {
6818: PetscInt i;
6821: if (!*mat) return(0);
6822: if (n < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Trying to destroy negative number of matrices %D",n);
6825: for (i=0; i<n; i++) {
6826: MatDestroy(&(*mat)[i]);
6827: }
6829: /* memory is allocated even if n = 0 */
6830: PetscFree(*mat);
6831: return(0);
6832: }
6834: /*@C
6835: MatDestroySubMatrices - Destroys a set of matrices obtained with MatCreateSubMatrices().
6837: Collective on Mat
6839: Input Parameters:
6840: + n - the number of local matrices
6841: - mat - the matrices (note that this is a pointer to the array of matrices, just to match the calling
6842: sequence of MatCreateSubMatrices())
6844: Level: advanced
6846: Notes: Frees not only the matrices, but also the array that contains the matrices
6847: In Fortran will not free the array.
6849: .seealso: MatCreateSubMatrices()
6850: @*/
6851: PetscErrorCode MatDestroySubMatrices(PetscInt n,Mat *mat[])
6852: {
6854: Mat mat0;
6857: if (!*mat) return(0);
6858: /* mat[] is an array of length n+1, see MatCreateSubMatrices_xxx() */
6859: if (n < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Trying to destroy negative number of matrices %D",n);
6862: mat0 = (*mat)[0];
6863: if (mat0 && mat0->ops->destroysubmatrices) {
6864: (mat0->ops->destroysubmatrices)(n,mat);
6865: } else {
6866: MatDestroyMatrices(n,mat);
6867: }
6868: return(0);
6869: }
6871: /*@C
6872: MatGetSeqNonzeroStructure - Extracts the sequential nonzero structure from a matrix.
6874: Collective on Mat
6876: Input Parameters:
6877: . mat - the matrix
6879: Output Parameter:
6880: . matstruct - the sequential matrix with the nonzero structure of mat
6882: Level: intermediate
6884: .seealso: MatDestroySeqNonzeroStructure(), MatCreateSubMatrices(), MatDestroyMatrices()
6885: @*/
6886: PetscErrorCode MatGetSeqNonzeroStructure(Mat mat,Mat *matstruct)
6887: {
6895: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
6896: MatCheckPreallocated(mat,1);
6898: if (!mat->ops->getseqnonzerostructure) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Not for matrix type %s\n",((PetscObject)mat)->type_name);
6899: PetscLogEventBegin(MAT_GetSeqNonzeroStructure,mat,0,0,0);
6900: (*mat->ops->getseqnonzerostructure)(mat,matstruct);
6901: PetscLogEventEnd(MAT_GetSeqNonzeroStructure,mat,0,0,0);
6902: return(0);
6903: }
6905: /*@C
6906: MatDestroySeqNonzeroStructure - Destroys matrix obtained with MatGetSeqNonzeroStructure().
6908: Collective on Mat
6910: Input Parameters:
6911: . mat - the matrix (note that this is a pointer to the array of matrices, just to match the calling
6912: sequence of MatGetSequentialNonzeroStructure())
6914: Level: advanced
6916: Notes: Frees not only the matrices, but also the array that contains the matrices
6918: .seealso: MatGetSeqNonzeroStructure()
6919: @*/
6920: PetscErrorCode MatDestroySeqNonzeroStructure(Mat *mat)
6921: {
6926: MatDestroy(mat);
6927: return(0);
6928: }
6930: /*@
6931: MatIncreaseOverlap - Given a set of submatrices indicated by index sets,
6932: replaces the index sets by larger ones that represent submatrices with
6933: additional overlap.
6935: Collective on Mat
6937: Input Parameters:
6938: + mat - the matrix
6939: . n - the number of index sets
6940: . is - the array of index sets (these index sets will changed during the call)
6941: - ov - the additional overlap requested
6943: Options Database:
6944: . -mat_increase_overlap_scalable - use a scalable algorithm to compute the overlap (supported by MPIAIJ matrix)
6946: Level: developer
6948: Concepts: overlap
6949: Concepts: ASM^computing overlap
6951: .seealso: MatCreateSubMatrices()
6952: @*/
6953: PetscErrorCode MatIncreaseOverlap(Mat mat,PetscInt n,IS is[],PetscInt ov)
6954: {
6960: if (n < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must have one or more domains, you have %D",n);
6961: if (n) {
6964: }
6965: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
6966: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
6967: MatCheckPreallocated(mat,1);
6969: if (!ov) return(0);
6970: if (!mat->ops->increaseoverlap) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
6971: PetscLogEventBegin(MAT_IncreaseOverlap,mat,0,0,0);
6972: (*mat->ops->increaseoverlap)(mat,n,is,ov);
6973: PetscLogEventEnd(MAT_IncreaseOverlap,mat,0,0,0);
6974: return(0);
6975: }
6978: PetscErrorCode MatIncreaseOverlapSplit_Single(Mat,IS*,PetscInt);
6980: /*@
6981: MatIncreaseOverlapSplit - Given a set of submatrices indicated by index sets across
6982: a sub communicator, replaces the index sets by larger ones that represent submatrices with
6983: additional overlap.
6985: Collective on Mat
6987: Input Parameters:
6988: + mat - the matrix
6989: . n - the number of index sets
6990: . is - the array of index sets (these index sets will changed during the call)
6991: - ov - the additional overlap requested
6993: Options Database:
6994: . -mat_increase_overlap_scalable - use a scalable algorithm to compute the overlap (supported by MPIAIJ matrix)
6996: Level: developer
6998: Concepts: overlap
6999: Concepts: ASM^computing overlap
7001: .seealso: MatCreateSubMatrices()
7002: @*/
7003: PetscErrorCode MatIncreaseOverlapSplit(Mat mat,PetscInt n,IS is[],PetscInt ov)
7004: {
7005: PetscInt i;
7011: if (n < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must have one or more domains, you have %D",n);
7012: if (n) {
7015: }
7016: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
7017: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
7018: MatCheckPreallocated(mat,1);
7019: if (!ov) return(0);
7020: PetscLogEventBegin(MAT_IncreaseOverlap,mat,0,0,0);
7021: for(i=0; i<n; i++){
7022: MatIncreaseOverlapSplit_Single(mat,&is[i],ov);
7023: }
7024: PetscLogEventEnd(MAT_IncreaseOverlap,mat,0,0,0);
7025: return(0);
7026: }
7031: /*@
7032: MatGetBlockSize - Returns the matrix block size.
7034: Not Collective
7036: Input Parameter:
7037: . mat - the matrix
7039: Output Parameter:
7040: . bs - block size
7042: Notes:
7043: Block row formats are MATSEQBAIJ, MATMPIBAIJ, MATSEQSBAIJ, MATMPISBAIJ. These formats ALWAYS have square block storage in the matrix.
7045: If the block size has not been set yet this routine returns 1.
7047: Level: intermediate
7049: Concepts: matrices^block size
7051: .seealso: MatCreateSeqBAIJ(), MatCreateBAIJ(), MatGetBlockSizes()
7052: @*/
7053: PetscErrorCode MatGetBlockSize(Mat mat,PetscInt *bs)
7054: {
7058: *bs = PetscAbs(mat->rmap->bs);
7059: return(0);
7060: }
7062: /*@
7063: MatGetBlockSizes - Returns the matrix block row and column sizes.
7065: Not Collective
7067: Input Parameter:
7068: . mat - the matrix
7070: Output Parameter:
7071: . rbs - row block size
7072: . cbs - column block size
7074: Notes:
7075: Block row formats are MATSEQBAIJ, MATMPIBAIJ, MATSEQSBAIJ, MATMPISBAIJ. These formats ALWAYS have square block storage in the matrix.
7076: If you pass a different block size for the columns than the rows, the row block size determines the square block storage.
7078: If a block size has not been set yet this routine returns 1.
7080: Level: intermediate
7082: Concepts: matrices^block size
7084: .seealso: MatCreateSeqBAIJ(), MatCreateBAIJ(), MatGetBlockSize(), MatSetBlockSize(), MatSetBlockSizes()
7085: @*/
7086: PetscErrorCode MatGetBlockSizes(Mat mat,PetscInt *rbs, PetscInt *cbs)
7087: {
7092: if (rbs) *rbs = PetscAbs(mat->rmap->bs);
7093: if (cbs) *cbs = PetscAbs(mat->cmap->bs);
7094: return(0);
7095: }
7097: /*@
7098: MatSetBlockSize - Sets the matrix block size.
7100: Logically Collective on Mat
7102: Input Parameters:
7103: + mat - the matrix
7104: - bs - block size
7106: Notes:
7107: Block row formats are MATSEQBAIJ, MATMPIBAIJ, MATSEQSBAIJ, MATMPISBAIJ. These formats ALWAYS have square block storage in the matrix.
7108: This must be called before MatSetUp() or MatXXXSetPreallocation() (or will default to 1) and the block size cannot be changed later.
7110: For MATMPIAIJ and MATSEQAIJ matrix formats, this function can be called at a later stage, provided that the specified block size
7111: is compatible with the matrix local sizes.
7113: Level: intermediate
7115: Concepts: matrices^block size
7117: .seealso: MatCreateSeqBAIJ(), MatCreateBAIJ(), MatGetBlockSize(), MatSetBlockSizes(), MatGetBlockSizes()
7118: @*/
7119: PetscErrorCode MatSetBlockSize(Mat mat,PetscInt bs)
7120: {
7126: MatSetBlockSizes(mat,bs,bs);
7127: return(0);
7128: }
7130: /*@
7131: MatSetBlockSizes - Sets the matrix block row and column sizes.
7133: Logically Collective on Mat
7135: Input Parameters:
7136: + mat - the matrix
7137: - rbs - row block size
7138: - cbs - column block size
7140: Notes:
7141: Block row formats are MATSEQBAIJ, MATMPIBAIJ, MATSEQSBAIJ, MATMPISBAIJ. These formats ALWAYS have square block storage in the matrix.
7142: If you pass a different block size for the columns than the rows, the row block size determines the square block storage.
7143: This must be called before MatSetUp() or MatXXXSetPreallocation() (or will default to 1) and the block size cannot be changed later
7145: For MATMPIAIJ and MATSEQAIJ matrix formats, this function can be called at a later stage, provided that the specified block sizes
7146: are compatible with the matrix local sizes.
7148: The row and column block size determine the blocksize of the "row" and "column" vectors returned by MatCreateVecs().
7150: Level: intermediate
7152: Concepts: matrices^block size
7154: .seealso: MatCreateSeqBAIJ(), MatCreateBAIJ(), MatGetBlockSize(), MatSetBlockSize(), MatGetBlockSizes()
7155: @*/
7156: PetscErrorCode MatSetBlockSizes(Mat mat,PetscInt rbs,PetscInt cbs)
7157: {
7164: if (mat->ops->setblocksizes) {
7165: (*mat->ops->setblocksizes)(mat,rbs,cbs);
7166: }
7167: if (mat->rmap->refcnt) {
7168: ISLocalToGlobalMapping l2g = NULL;
7169: PetscLayout nmap = NULL;
7171: PetscLayoutDuplicate(mat->rmap,&nmap);
7172: if (mat->rmap->mapping) {
7173: ISLocalToGlobalMappingDuplicate(mat->rmap->mapping,&l2g);
7174: }
7175: PetscLayoutDestroy(&mat->rmap);
7176: mat->rmap = nmap;
7177: mat->rmap->mapping = l2g;
7178: }
7179: if (mat->cmap->refcnt) {
7180: ISLocalToGlobalMapping l2g = NULL;
7181: PetscLayout nmap = NULL;
7183: PetscLayoutDuplicate(mat->cmap,&nmap);
7184: if (mat->cmap->mapping) {
7185: ISLocalToGlobalMappingDuplicate(mat->cmap->mapping,&l2g);
7186: }
7187: PetscLayoutDestroy(&mat->cmap);
7188: mat->cmap = nmap;
7189: mat->cmap->mapping = l2g;
7190: }
7191: PetscLayoutSetBlockSize(mat->rmap,rbs);
7192: PetscLayoutSetBlockSize(mat->cmap,cbs);
7193: return(0);
7194: }
7196: /*@
7197: MatSetBlockSizesFromMats - Sets the matrix block row and column sizes to match a pair of matrices
7199: Logically Collective on Mat
7201: Input Parameters:
7202: + mat - the matrix
7203: . fromRow - matrix from which to copy row block size
7204: - fromCol - matrix from which to copy column block size (can be same as fromRow)
7206: Level: developer
7208: Concepts: matrices^block size
7210: .seealso: MatCreateSeqBAIJ(), MatCreateBAIJ(), MatGetBlockSize(), MatSetBlockSizes()
7211: @*/
7212: PetscErrorCode MatSetBlockSizesFromMats(Mat mat,Mat fromRow,Mat fromCol)
7213: {
7220: if (fromRow->rmap->bs > 0) {PetscLayoutSetBlockSize(mat->rmap,fromRow->rmap->bs);}
7221: if (fromCol->cmap->bs > 0) {PetscLayoutSetBlockSize(mat->cmap,fromCol->cmap->bs);}
7222: return(0);
7223: }
7225: /*@
7226: MatResidual - Default routine to calculate the residual.
7228: Collective on Mat and Vec
7230: Input Parameters:
7231: + mat - the matrix
7232: . b - the right-hand-side
7233: - x - the approximate solution
7235: Output Parameter:
7236: . r - location to store the residual
7238: Level: developer
7240: .keywords: MG, default, multigrid, residual
7242: .seealso: PCMGSetResidual()
7243: @*/
7244: PetscErrorCode MatResidual(Mat mat,Vec b,Vec x,Vec r)
7245: {
7254: MatCheckPreallocated(mat,1);
7255: PetscLogEventBegin(MAT_Residual,mat,0,0,0);
7256: if (!mat->ops->residual) {
7257: MatMult(mat,x,r);
7258: VecAYPX(r,-1.0,b);
7259: } else {
7260: (*mat->ops->residual)(mat,b,x,r);
7261: }
7262: PetscLogEventEnd(MAT_Residual,mat,0,0,0);
7263: return(0);
7264: }
7266: /*@C
7267: MatGetRowIJ - Returns the compressed row storage i and j indices for sequential matrices.
7269: Collective on Mat
7271: Input Parameters:
7272: + mat - the matrix
7273: . shift - 0 or 1 indicating we want the indices starting at 0 or 1
7274: . symmetric - PETSC_TRUE or PETSC_FALSE indicating the matrix data structure should be symmetrized
7275: - inodecompressed - PETSC_TRUE or PETSC_FALSE indicating if the nonzero structure of the
7276: inodes or the nonzero elements is wanted. For BAIJ matrices the compressed version is
7277: always used.
7279: Output Parameters:
7280: + n - number of rows in the (possibly compressed) matrix
7281: . ia - the row pointers [of length n+1]
7282: . ja - the column indices
7283: - done - indicates if the routine actually worked and returned appropriate ia[] and ja[] arrays; callers
7284: are responsible for handling the case when done == PETSC_FALSE and ia and ja are not set
7286: Level: developer
7288: Notes:
7289: You CANNOT change any of the ia[] or ja[] values.
7291: Use MatRestoreRowIJ() when you are finished accessing the ia[] and ja[] values.
7293: Fortran Notes:
7294: In Fortran use
7295: $
7296: $ PetscInt ia(1), ja(1)
7297: $ PetscOffset iia, jja
7298: $ call MatGetRowIJ(mat,shift,symmetric,inodecompressed,n,ia,iia,ja,jja,done,ierr)
7299: $ ! Access the ith and jth entries via ia(iia + i) and ja(jja + j)
7301: or
7302: $
7303: $ PetscInt, pointer :: ia(:),ja(:)
7304: $ call MatGetRowIJF90(mat,shift,symmetric,inodecompressed,n,ia,ja,done,ierr)
7305: $ ! Access the ith and jth entries via ia(i) and ja(j)
7307: .seealso: MatGetColumnIJ(), MatRestoreRowIJ(), MatSeqAIJGetArray()
7308: @*/
7309: PetscErrorCode MatGetRowIJ(Mat mat,PetscInt shift,PetscBool symmetric,PetscBool inodecompressed,PetscInt *n,const PetscInt *ia[],const PetscInt *ja[],PetscBool *done)
7310: {
7320: MatCheckPreallocated(mat,1);
7321: if (!mat->ops->getrowij) *done = PETSC_FALSE;
7322: else {
7323: *done = PETSC_TRUE;
7324: PetscLogEventBegin(MAT_GetRowIJ,mat,0,0,0);
7325: (*mat->ops->getrowij)(mat,shift,symmetric,inodecompressed,n,ia,ja,done);
7326: PetscLogEventEnd(MAT_GetRowIJ,mat,0,0,0);
7327: }
7328: return(0);
7329: }
7331: /*@C
7332: MatGetColumnIJ - Returns the compressed column storage i and j indices for sequential matrices.
7334: Collective on Mat
7336: Input Parameters:
7337: + mat - the matrix
7338: . shift - 1 or zero indicating we want the indices starting at 0 or 1
7339: . symmetric - PETSC_TRUE or PETSC_FALSE indicating the matrix data structure should be
7340: symmetrized
7341: . inodecompressed - PETSC_TRUE or PETSC_FALSE indicating if the nonzero structure of the
7342: inodes or the nonzero elements is wanted. For BAIJ matrices the compressed version is
7343: always used.
7344: . n - number of columns in the (possibly compressed) matrix
7345: . ia - the column pointers
7346: - ja - the row indices
7348: Output Parameters:
7349: . done - PETSC_TRUE or PETSC_FALSE, indicating whether the values have been returned
7351: Note:
7352: This routine zeros out n, ia, and ja. This is to prevent accidental
7353: us of the array after it has been restored. If you pass NULL, it will
7354: not zero the pointers. Use of ia or ja after MatRestoreColumnIJ() is invalid.
7356: Level: developer
7358: .seealso: MatGetRowIJ(), MatRestoreColumnIJ()
7359: @*/
7360: PetscErrorCode MatGetColumnIJ(Mat mat,PetscInt shift,PetscBool symmetric,PetscBool inodecompressed,PetscInt *n,const PetscInt *ia[],const PetscInt *ja[],PetscBool *done)
7361: {
7371: MatCheckPreallocated(mat,1);
7372: if (!mat->ops->getcolumnij) *done = PETSC_FALSE;
7373: else {
7374: *done = PETSC_TRUE;
7375: (*mat->ops->getcolumnij)(mat,shift,symmetric,inodecompressed,n,ia,ja,done);
7376: }
7377: return(0);
7378: }
7380: /*@C
7381: MatRestoreRowIJ - Call after you are completed with the ia,ja indices obtained with
7382: MatGetRowIJ().
7384: Collective on Mat
7386: Input Parameters:
7387: + mat - the matrix
7388: . shift - 1 or zero indicating we want the indices starting at 0 or 1
7389: . symmetric - PETSC_TRUE or PETSC_FALSE indicating the matrix data structure should be
7390: symmetrized
7391: . inodecompressed - PETSC_TRUE or PETSC_FALSE indicating if the nonzero structure of the
7392: inodes or the nonzero elements is wanted. For BAIJ matrices the compressed version is
7393: always used.
7394: . n - size of (possibly compressed) matrix
7395: . ia - the row pointers
7396: - ja - the column indices
7398: Output Parameters:
7399: . done - PETSC_TRUE or PETSC_FALSE indicated that the values have been returned
7401: Note:
7402: This routine zeros out n, ia, and ja. This is to prevent accidental
7403: us of the array after it has been restored. If you pass NULL, it will
7404: not zero the pointers. Use of ia or ja after MatRestoreRowIJ() is invalid.
7406: Level: developer
7408: .seealso: MatGetRowIJ(), MatRestoreColumnIJ()
7409: @*/
7410: PetscErrorCode MatRestoreRowIJ(Mat mat,PetscInt shift,PetscBool symmetric,PetscBool inodecompressed,PetscInt *n,const PetscInt *ia[],const PetscInt *ja[],PetscBool *done)
7411: {
7420: MatCheckPreallocated(mat,1);
7422: if (!mat->ops->restorerowij) *done = PETSC_FALSE;
7423: else {
7424: *done = PETSC_TRUE;
7425: (*mat->ops->restorerowij)(mat,shift,symmetric,inodecompressed,n,ia,ja,done);
7426: if (n) *n = 0;
7427: if (ia) *ia = NULL;
7428: if (ja) *ja = NULL;
7429: }
7430: return(0);
7431: }
7433: /*@C
7434: MatRestoreColumnIJ - Call after you are completed with the ia,ja indices obtained with
7435: MatGetColumnIJ().
7437: Collective on Mat
7439: Input Parameters:
7440: + mat - the matrix
7441: . shift - 1 or zero indicating we want the indices starting at 0 or 1
7442: - symmetric - PETSC_TRUE or PETSC_FALSE indicating the matrix data structure should be
7443: symmetrized
7444: - inodecompressed - PETSC_TRUE or PETSC_FALSE indicating if the nonzero structure of the
7445: inodes or the nonzero elements is wanted. For BAIJ matrices the compressed version is
7446: always used.
7448: Output Parameters:
7449: + n - size of (possibly compressed) matrix
7450: . ia - the column pointers
7451: . ja - the row indices
7452: - done - PETSC_TRUE or PETSC_FALSE indicated that the values have been returned
7454: Level: developer
7456: .seealso: MatGetColumnIJ(), MatRestoreRowIJ()
7457: @*/
7458: PetscErrorCode MatRestoreColumnIJ(Mat mat,PetscInt shift,PetscBool symmetric,PetscBool inodecompressed,PetscInt *n,const PetscInt *ia[],const PetscInt *ja[],PetscBool *done)
7459: {
7468: MatCheckPreallocated(mat,1);
7470: if (!mat->ops->restorecolumnij) *done = PETSC_FALSE;
7471: else {
7472: *done = PETSC_TRUE;
7473: (*mat->ops->restorecolumnij)(mat,shift,symmetric,inodecompressed,n,ia,ja,done);
7474: if (n) *n = 0;
7475: if (ia) *ia = NULL;
7476: if (ja) *ja = NULL;
7477: }
7478: return(0);
7479: }
7481: /*@C
7482: MatColoringPatch -Used inside matrix coloring routines that
7483: use MatGetRowIJ() and/or MatGetColumnIJ().
7485: Collective on Mat
7487: Input Parameters:
7488: + mat - the matrix
7489: . ncolors - max color value
7490: . n - number of entries in colorarray
7491: - colorarray - array indicating color for each column
7493: Output Parameters:
7494: . iscoloring - coloring generated using colorarray information
7496: Level: developer
7498: .seealso: MatGetRowIJ(), MatGetColumnIJ()
7500: @*/
7501: PetscErrorCode MatColoringPatch(Mat mat,PetscInt ncolors,PetscInt n,ISColoringValue colorarray[],ISColoring *iscoloring)
7502: {
7510: MatCheckPreallocated(mat,1);
7512: if (!mat->ops->coloringpatch) {
7513: ISColoringCreate(PetscObjectComm((PetscObject)mat),ncolors,n,colorarray,PETSC_OWN_POINTER,iscoloring);
7514: } else {
7515: (*mat->ops->coloringpatch)(mat,ncolors,n,colorarray,iscoloring);
7516: }
7517: return(0);
7518: }
7521: /*@
7522: MatSetUnfactored - Resets a factored matrix to be treated as unfactored.
7524: Logically Collective on Mat
7526: Input Parameter:
7527: . mat - the factored matrix to be reset
7529: Notes:
7530: This routine should be used only with factored matrices formed by in-place
7531: factorization via ILU(0) (or by in-place LU factorization for the MATSEQDENSE
7532: format). This option can save memory, for example, when solving nonlinear
7533: systems with a matrix-free Newton-Krylov method and a matrix-based, in-place
7534: ILU(0) preconditioner.
7536: Note that one can specify in-place ILU(0) factorization by calling
7537: .vb
7538: PCType(pc,PCILU);
7539: PCFactorSeUseInPlace(pc);
7540: .ve
7541: or by using the options -pc_type ilu -pc_factor_in_place
7543: In-place factorization ILU(0) can also be used as a local
7544: solver for the blocks within the block Jacobi or additive Schwarz
7545: methods (runtime option: -sub_pc_factor_in_place). See Users-Manual: ch_pc
7546: for details on setting local solver options.
7548: Most users should employ the simplified KSP interface for linear solvers
7549: instead of working directly with matrix algebra routines such as this.
7550: See, e.g., KSPCreate().
7552: Level: developer
7554: .seealso: PCFactorSetUseInPlace(), PCFactorGetUseInPlace()
7556: Concepts: matrices^unfactored
7558: @*/
7559: PetscErrorCode MatSetUnfactored(Mat mat)
7560: {
7566: MatCheckPreallocated(mat,1);
7567: mat->factortype = MAT_FACTOR_NONE;
7568: if (!mat->ops->setunfactored) return(0);
7569: (*mat->ops->setunfactored)(mat);
7570: return(0);
7571: }
7573: /*MC
7574: MatDenseGetArrayF90 - Accesses a matrix array from Fortran90.
7576: Synopsis:
7577: MatDenseGetArrayF90(Mat x,{Scalar, pointer :: xx_v(:,:)},integer ierr)
7579: Not collective
7581: Input Parameter:
7582: . x - matrix
7584: Output Parameters:
7585: + xx_v - the Fortran90 pointer to the array
7586: - ierr - error code
7588: Example of Usage:
7589: .vb
7590: PetscScalar, pointer xx_v(:,:)
7591: ....
7592: call MatDenseGetArrayF90(x,xx_v,ierr)
7593: a = xx_v(3)
7594: call MatDenseRestoreArrayF90(x,xx_v,ierr)
7595: .ve
7597: Level: advanced
7599: .seealso: MatDenseRestoreArrayF90(), MatDenseGetArray(), MatDenseRestoreArray(), MatSeqAIJGetArrayF90()
7601: Concepts: matrices^accessing array
7603: M*/
7605: /*MC
7606: MatDenseRestoreArrayF90 - Restores a matrix array that has been
7607: accessed with MatDenseGetArrayF90().
7609: Synopsis:
7610: MatDenseRestoreArrayF90(Mat x,{Scalar, pointer :: xx_v(:,:)},integer ierr)
7612: Not collective
7614: Input Parameters:
7615: + x - matrix
7616: - xx_v - the Fortran90 pointer to the array
7618: Output Parameter:
7619: . ierr - error code
7621: Example of Usage:
7622: .vb
7623: PetscScalar, pointer xx_v(:,:)
7624: ....
7625: call MatDenseGetArrayF90(x,xx_v,ierr)
7626: a = xx_v(3)
7627: call MatDenseRestoreArrayF90(x,xx_v,ierr)
7628: .ve
7630: Level: advanced
7632: .seealso: MatDenseGetArrayF90(), MatDenseGetArray(), MatDenseRestoreArray(), MatSeqAIJRestoreArrayF90()
7634: M*/
7637: /*MC
7638: MatSeqAIJGetArrayF90 - Accesses a matrix array from Fortran90.
7640: Synopsis:
7641: MatSeqAIJGetArrayF90(Mat x,{Scalar, pointer :: xx_v(:)},integer ierr)
7643: Not collective
7645: Input Parameter:
7646: . x - matrix
7648: Output Parameters:
7649: + xx_v - the Fortran90 pointer to the array
7650: - ierr - error code
7652: Example of Usage:
7653: .vb
7654: PetscScalar, pointer xx_v(:)
7655: ....
7656: call MatSeqAIJGetArrayF90(x,xx_v,ierr)
7657: a = xx_v(3)
7658: call MatSeqAIJRestoreArrayF90(x,xx_v,ierr)
7659: .ve
7661: Level: advanced
7663: .seealso: MatSeqAIJRestoreArrayF90(), MatSeqAIJGetArray(), MatSeqAIJRestoreArray(), MatDenseGetArrayF90()
7665: Concepts: matrices^accessing array
7667: M*/
7669: /*MC
7670: MatSeqAIJRestoreArrayF90 - Restores a matrix array that has been
7671: accessed with MatSeqAIJGetArrayF90().
7673: Synopsis:
7674: MatSeqAIJRestoreArrayF90(Mat x,{Scalar, pointer :: xx_v(:)},integer ierr)
7676: Not collective
7678: Input Parameters:
7679: + x - matrix
7680: - xx_v - the Fortran90 pointer to the array
7682: Output Parameter:
7683: . ierr - error code
7685: Example of Usage:
7686: .vb
7687: PetscScalar, pointer xx_v(:)
7688: ....
7689: call MatSeqAIJGetArrayF90(x,xx_v,ierr)
7690: a = xx_v(3)
7691: call MatSeqAIJRestoreArrayF90(x,xx_v,ierr)
7692: .ve
7694: Level: advanced
7696: .seealso: MatSeqAIJGetArrayF90(), MatSeqAIJGetArray(), MatSeqAIJRestoreArray(), MatDenseRestoreArrayF90()
7698: M*/
7701: /*@
7702: MatCreateSubMatrix - Gets a single submatrix on the same number of processors
7703: as the original matrix.
7705: Collective on Mat
7707: Input Parameters:
7708: + mat - the original matrix
7709: . isrow - parallel IS containing the rows this processor should obtain
7710: . iscol - parallel IS containing all columns you wish to keep. Each process should list the columns that will be in IT's "diagonal part" in the new matrix.
7711: - cll - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
7713: Output Parameter:
7714: . newmat - the new submatrix, of the same type as the old
7716: Level: advanced
7718: Notes:
7719: The submatrix will be able to be multiplied with vectors using the same layout as iscol.
7721: Some matrix types place restrictions on the row and column indices, such
7722: as that they be sorted or that they be equal to each other.
7724: The index sets may not have duplicate entries.
7726: The first time this is called you should use a cll of MAT_INITIAL_MATRIX,
7727: the MatCreateSubMatrix() routine will create the newmat for you. Any additional calls
7728: to this routine with a mat of the same nonzero structure and with a call of MAT_REUSE_MATRIX
7729: will reuse the matrix generated the first time. You should call MatDestroy() on newmat when
7730: you are finished using it.
7732: The communicator of the newly obtained matrix is ALWAYS the same as the communicator of
7733: the input matrix.
7735: If iscol is NULL then all columns are obtained (not supported in Fortran).
7737: Example usage:
7738: Consider the following 8x8 matrix with 34 non-zero values, that is
7739: assembled across 3 processors. Let's assume that proc0 owns 3 rows,
7740: proc1 owns 3 rows, proc2 owns 2 rows. This division can be shown
7741: as follows:
7743: .vb
7744: 1 2 0 | 0 3 0 | 0 4
7745: Proc0 0 5 6 | 7 0 0 | 8 0
7746: 9 0 10 | 11 0 0 | 12 0
7747: -------------------------------------
7748: 13 0 14 | 15 16 17 | 0 0
7749: Proc1 0 18 0 | 19 20 21 | 0 0
7750: 0 0 0 | 22 23 0 | 24 0
7751: -------------------------------------
7752: Proc2 25 26 27 | 0 0 28 | 29 0
7753: 30 0 0 | 31 32 33 | 0 34
7754: .ve
7756: Suppose isrow = [0 1 | 4 | 6 7] and iscol = [1 2 | 3 4 5 | 6]. The resulting submatrix is
7758: .vb
7759: 2 0 | 0 3 0 | 0
7760: Proc0 5 6 | 7 0 0 | 8
7761: -------------------------------
7762: Proc1 18 0 | 19 20 21 | 0
7763: -------------------------------
7764: Proc2 26 27 | 0 0 28 | 29
7765: 0 0 | 31 32 33 | 0
7766: .ve
7769: Concepts: matrices^submatrices
7771: .seealso: MatCreateSubMatrices()
7772: @*/
7773: PetscErrorCode MatCreateSubMatrix(Mat mat,IS isrow,IS iscol,MatReuse cll,Mat *newmat)
7774: {
7776: PetscMPIInt size;
7777: Mat *local;
7778: IS iscoltmp;
7787: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
7788: if (cll == MAT_IGNORE_MATRIX) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Cannot use MAT_IGNORE_MATRIX");
7790: MatCheckPreallocated(mat,1);
7791: MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);
7793: if (!iscol || isrow == iscol) {
7794: PetscBool stride;
7795: PetscMPIInt grabentirematrix = 0,grab;
7796: PetscObjectTypeCompare((PetscObject)isrow,ISSTRIDE,&stride);
7797: if (stride) {
7798: PetscInt first,step,n,rstart,rend;
7799: ISStrideGetInfo(isrow,&first,&step);
7800: if (step == 1) {
7801: MatGetOwnershipRange(mat,&rstart,&rend);
7802: if (rstart == first) {
7803: ISGetLocalSize(isrow,&n);
7804: if (n == rend-rstart) {
7805: grabentirematrix = 1;
7806: }
7807: }
7808: }
7809: }
7810: MPIU_Allreduce(&grabentirematrix,&grab,1,MPI_INT,MPI_MIN,PetscObjectComm((PetscObject)mat));
7811: if (grab) {
7812: PetscInfo(mat,"Getting entire matrix as submatrix\n");
7813: if (cll == MAT_INITIAL_MATRIX) {
7814: *newmat = mat;
7815: PetscObjectReference((PetscObject)mat);
7816: }
7817: return(0);
7818: }
7819: }
7821: if (!iscol) {
7822: ISCreateStride(PetscObjectComm((PetscObject)mat),mat->cmap->n,mat->cmap->rstart,1,&iscoltmp);
7823: } else {
7824: iscoltmp = iscol;
7825: }
7827: /* if original matrix is on just one processor then use submatrix generated */
7828: if (mat->ops->createsubmatrices && !mat->ops->createsubmatrix && size == 1 && cll == MAT_REUSE_MATRIX) {
7829: MatCreateSubMatrices(mat,1,&isrow,&iscoltmp,MAT_REUSE_MATRIX,&newmat);
7830: if (!iscol) {ISDestroy(&iscoltmp);}
7831: return(0);
7832: } else if (mat->ops->createsubmatrices && !mat->ops->createsubmatrix && size == 1) {
7833: MatCreateSubMatrices(mat,1,&isrow,&iscoltmp,MAT_INITIAL_MATRIX,&local);
7834: *newmat = *local;
7835: PetscFree(local);
7836: if (!iscol) {ISDestroy(&iscoltmp);}
7837: return(0);
7838: } else if (!mat->ops->createsubmatrix) {
7839: /* Create a new matrix type that implements the operation using the full matrix */
7840: PetscLogEventBegin(MAT_CreateSubMat,mat,0,0,0);
7841: switch (cll) {
7842: case MAT_INITIAL_MATRIX:
7843: MatCreateSubMatrixVirtual(mat,isrow,iscoltmp,newmat);
7844: break;
7845: case MAT_REUSE_MATRIX:
7846: MatSubMatrixVirtualUpdate(*newmat,mat,isrow,iscoltmp);
7847: break;
7848: default: SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_OUTOFRANGE,"Invalid MatReuse, must be either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX");
7849: }
7850: PetscLogEventEnd(MAT_CreateSubMat,mat,0,0,0);
7851: if (!iscol) {ISDestroy(&iscoltmp);}
7852: return(0);
7853: }
7855: if (!mat->ops->createsubmatrix) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
7856: PetscLogEventBegin(MAT_CreateSubMat,mat,0,0,0);
7857: (*mat->ops->createsubmatrix)(mat,isrow,iscoltmp,cll,newmat);
7858: PetscLogEventEnd(MAT_CreateSubMat,mat,0,0,0);
7860: /* Propagate symmetry information for diagonal blocks */
7861: if (isrow == iscoltmp) {
7862: if (mat->symmetric_set && mat->symmetric) {
7863: MatSetOption(*newmat,MAT_SYMMETRIC,PETSC_TRUE);
7864: }
7865: if (mat->structurally_symmetric_set && mat->structurally_symmetric) {
7866: MatSetOption(*newmat,MAT_STRUCTURALLY_SYMMETRIC,PETSC_TRUE);
7867: }
7868: if (mat->hermitian_set && mat->hermitian) {
7869: MatSetOption(*newmat,MAT_HERMITIAN,PETSC_TRUE);
7870: }
7871: if (mat->spd_set && mat->spd) {
7872: MatSetOption(*newmat,MAT_SPD,PETSC_TRUE);
7873: }
7874: }
7876: if (!iscol) {ISDestroy(&iscoltmp);}
7877: if (*newmat && cll == MAT_INITIAL_MATRIX) {PetscObjectStateIncrease((PetscObject)*newmat);}
7878: return(0);
7879: }
7881: /*@
7882: MatStashSetInitialSize - sets the sizes of the matrix stash, that is
7883: used during the assembly process to store values that belong to
7884: other processors.
7886: Not Collective
7888: Input Parameters:
7889: + mat - the matrix
7890: . size - the initial size of the stash.
7891: - bsize - the initial size of the block-stash(if used).
7893: Options Database Keys:
7894: + -matstash_initial_size <size> or <size0,size1,...sizep-1>
7895: - -matstash_block_initial_size <bsize> or <bsize0,bsize1,...bsizep-1>
7897: Level: intermediate
7899: Notes:
7900: The block-stash is used for values set with MatSetValuesBlocked() while
7901: the stash is used for values set with MatSetValues()
7903: Run with the option -info and look for output of the form
7904: MatAssemblyBegin_MPIXXX:Stash has MM entries, uses nn mallocs.
7905: to determine the appropriate value, MM, to use for size and
7906: MatAssemblyBegin_MPIXXX:Block-Stash has BMM entries, uses nn mallocs.
7907: to determine the value, BMM to use for bsize
7909: Concepts: stash^setting matrix size
7910: Concepts: matrices^stash
7912: .seealso: MatAssemblyBegin(), MatAssemblyEnd(), Mat, MatStashGetInfo()
7914: @*/
7915: PetscErrorCode MatStashSetInitialSize(Mat mat,PetscInt size, PetscInt bsize)
7916: {
7922: MatStashSetInitialSize_Private(&mat->stash,size);
7923: MatStashSetInitialSize_Private(&mat->bstash,bsize);
7924: return(0);
7925: }
7927: /*@
7928: MatInterpolateAdd - w = y + A*x or A'*x depending on the shape of
7929: the matrix
7931: Neighbor-wise Collective on Mat
7933: Input Parameters:
7934: + mat - the matrix
7935: . x,y - the vectors
7936: - w - where the result is stored
7938: Level: intermediate
7940: Notes:
7941: w may be the same vector as y.
7943: This allows one to use either the restriction or interpolation (its transpose)
7944: matrix to do the interpolation
7946: Concepts: interpolation
7948: .seealso: MatMultAdd(), MatMultTransposeAdd(), MatRestrict()
7950: @*/
7951: PetscErrorCode MatInterpolateAdd(Mat A,Vec x,Vec y,Vec w)
7952: {
7954: PetscInt M,N,Ny;
7962: MatCheckPreallocated(A,1);
7963: MatGetSize(A,&M,&N);
7964: VecGetSize(y,&Ny);
7965: if (M == Ny) {
7966: MatMultAdd(A,x,y,w);
7967: } else {
7968: MatMultTransposeAdd(A,x,y,w);
7969: }
7970: return(0);
7971: }
7973: /*@
7974: MatInterpolate - y = A*x or A'*x depending on the shape of
7975: the matrix
7977: Neighbor-wise Collective on Mat
7979: Input Parameters:
7980: + mat - the matrix
7981: - x,y - the vectors
7983: Level: intermediate
7985: Notes:
7986: This allows one to use either the restriction or interpolation (its transpose)
7987: matrix to do the interpolation
7989: Concepts: matrices^interpolation
7991: .seealso: MatMultAdd(), MatMultTransposeAdd(), MatRestrict()
7993: @*/
7994: PetscErrorCode MatInterpolate(Mat A,Vec x,Vec y)
7995: {
7997: PetscInt M,N,Ny;
8004: MatCheckPreallocated(A,1);
8005: MatGetSize(A,&M,&N);
8006: VecGetSize(y,&Ny);
8007: if (M == Ny) {
8008: MatMult(A,x,y);
8009: } else {
8010: MatMultTranspose(A,x,y);
8011: }
8012: return(0);
8013: }
8015: /*@
8016: MatRestrict - y = A*x or A'*x
8018: Neighbor-wise Collective on Mat
8020: Input Parameters:
8021: + mat - the matrix
8022: - x,y - the vectors
8024: Level: intermediate
8026: Notes:
8027: This allows one to use either the restriction or interpolation (its transpose)
8028: matrix to do the restriction
8030: Concepts: matrices^restriction
8032: .seealso: MatMultAdd(), MatMultTransposeAdd(), MatInterpolate()
8034: @*/
8035: PetscErrorCode MatRestrict(Mat A,Vec x,Vec y)
8036: {
8038: PetscInt M,N,Ny;
8045: MatCheckPreallocated(A,1);
8047: MatGetSize(A,&M,&N);
8048: VecGetSize(y,&Ny);
8049: if (M == Ny) {
8050: MatMult(A,x,y);
8051: } else {
8052: MatMultTranspose(A,x,y);
8053: }
8054: return(0);
8055: }
8057: /*@C
8058: MatGetNullSpace - retrieves the null space of a matrix.
8060: Logically Collective on Mat and MatNullSpace
8062: Input Parameters:
8063: + mat - the matrix
8064: - nullsp - the null space object
8066: Level: developer
8068: Concepts: null space^attaching to matrix
8070: .seealso: MatCreate(), MatNullSpaceCreate(), MatSetNearNullSpace(), MatSetNullSpace()
8071: @*/
8072: PetscErrorCode MatGetNullSpace(Mat mat, MatNullSpace *nullsp)
8073: {
8077: *nullsp = (mat->symmetric_set && mat->symmetric && !mat->nullsp) ? mat->transnullsp : mat->nullsp;
8078: return(0);
8079: }
8081: /*@C
8082: MatSetNullSpace - attaches a null space to a matrix.
8084: Logically Collective on Mat and MatNullSpace
8086: Input Parameters:
8087: + mat - the matrix
8088: - nullsp - the null space object
8090: Level: advanced
8092: Notes:
8093: This null space is used by the linear solvers. Overwrites any previous null space that may have been attached
8095: For inconsistent singular systems (linear systems where the right hand side is not in the range of the operator) you also likely should
8096: call MatSetTransposeNullSpace(). This allows the linear system to be solved in a least squares sense.
8098: You can remove the null space by calling this routine with an nullsp of NULL
8101: The fundamental theorem of linear algebra (Gilbert Strang, Introduction to Applied Mathematics, page 72) states that
8102: the domain of a matrix A (from R^n to R^m (m rows, n columns) R^n = the direct sum of the null space of A, n(A), + the range of A^T, R(A^T).
8103: Similarly R^m = direct sum n(A^T) + R(A). Hence the linear system A x = b has a solution only if b in R(A) (or correspondingly b is orthogonal to
8104: n(A^T)) and if x is a solution then x + alpha n(A) is a solution for any alpha. The minimum norm solution is orthogonal to n(A). For problems without a solution
8105: the solution that minimizes the norm of the residual (the least squares solution) can be obtained by solving A x = \hat{b} where \hat{b} is b orthogonalized to the n(A^T).
8107: Krylov solvers can produce the minimal norm solution to the least squares problem by utilizing MatNullSpaceRemove().
8109: If the matrix is known to be symmetric because it is an SBAIJ matrix or one as called MatSetOption(mat,MAT_SYMMETRIC or MAT_SYMMETRIC_ETERNAL,PETSC_TRUE); this
8110: routine also automatically calls MatSetTransposeNullSpace().
8112: Concepts: null space^attaching to matrix
8114: .seealso: MatCreate(), MatNullSpaceCreate(), MatSetNearNullSpace(), MatGetNullSpace(), MatSetTransposeNullSpace(), MatGetTransposeNullSpace(), MatNullSpaceRemove()
8115: @*/
8116: PetscErrorCode MatSetNullSpace(Mat mat,MatNullSpace nullsp)
8117: {
8123: if (nullsp) {PetscObjectReference((PetscObject)nullsp);}
8124: MatNullSpaceDestroy(&mat->nullsp);
8125: mat->nullsp = nullsp;
8126: if (mat->symmetric_set && mat->symmetric) {
8127: MatSetTransposeNullSpace(mat,nullsp);
8128: }
8129: return(0);
8130: }
8132: /*@
8133: MatGetTransposeNullSpace - retrieves the null space of the transpose of a matrix.
8135: Logically Collective on Mat and MatNullSpace
8137: Input Parameters:
8138: + mat - the matrix
8139: - nullsp - the null space object
8141: Level: developer
8143: Concepts: null space^attaching to matrix
8145: .seealso: MatCreate(), MatNullSpaceCreate(), MatSetNearNullSpace(), MatSetTransposeNullSpace(), MatSetNullSpace(), MatGetNullSpace()
8146: @*/
8147: PetscErrorCode MatGetTransposeNullSpace(Mat mat, MatNullSpace *nullsp)
8148: {
8153: *nullsp = (mat->symmetric_set && mat->symmetric && !mat->transnullsp) ? mat->nullsp : mat->transnullsp;
8154: return(0);
8155: }
8157: /*@
8158: MatSetTransposeNullSpace - attaches a null space to a matrix.
8160: Logically Collective on Mat and MatNullSpace
8162: Input Parameters:
8163: + mat - the matrix
8164: - nullsp - the null space object
8166: Level: advanced
8168: Notes:
8169: For inconsistent singular systems (linear systems where the right hand side is not in the range of the operator) this allows the linear system to be solved in a least squares sense.
8170: You must also call MatSetNullSpace()
8173: The fundamental theorem of linear algebra (Gilbert Strang, Introduction to Applied Mathematics, page 72) states that
8174: the domain of a matrix A (from R^n to R^m (m rows, n columns) R^n = the direct sum of the null space of A, n(A), + the range of A^T, R(A^T).
8175: Similarly R^m = direct sum n(A^T) + R(A). Hence the linear system A x = b has a solution only if b in R(A) (or correspondingly b is orthogonal to
8176: n(A^T)) and if x is a solution then x + alpha n(A) is a solution for any alpha. The minimum norm solution is orthogonal to n(A). For problems without a solution
8177: the solution that minimizes the norm of the residual (the least squares solution) can be obtained by solving A x = \hat{b} where \hat{b} is b orthogonalized to the n(A^T).
8179: Krylov solvers can produce the minimal norm solution to the least squares problem by utilizing MatNullSpaceRemove().
8181: Concepts: null space^attaching to matrix
8183: .seealso: MatCreate(), MatNullSpaceCreate(), MatSetNearNullSpace(), MatGetNullSpace(), MatSetNullSpace(), MatGetTransposeNullSpace(), MatNullSpaceRemove()
8184: @*/
8185: PetscErrorCode MatSetTransposeNullSpace(Mat mat,MatNullSpace nullsp)
8186: {
8193: MatCheckPreallocated(mat,1);
8194: PetscObjectReference((PetscObject)nullsp);
8195: MatNullSpaceDestroy(&mat->transnullsp);
8196: mat->transnullsp = nullsp;
8197: return(0);
8198: }
8200: /*@
8201: MatSetNearNullSpace - attaches a null space to a matrix, which is often the null space (rigid body modes) of the operator without boundary conditions
8202: This null space will be used to provide near null space vectors to a multigrid preconditioner built from this matrix.
8204: Logically Collective on Mat and MatNullSpace
8206: Input Parameters:
8207: + mat - the matrix
8208: - nullsp - the null space object
8210: Level: advanced
8212: Notes:
8213: Overwrites any previous near null space that may have been attached
8215: You can remove the null space by calling this routine with an nullsp of NULL
8217: Concepts: null space^attaching to matrix
8219: .seealso: MatCreate(), MatNullSpaceCreate(), MatSetNullSpace(), MatNullSpaceCreateRigidBody(), MatGetNearNullSpace()
8220: @*/
8221: PetscErrorCode MatSetNearNullSpace(Mat mat,MatNullSpace nullsp)
8222: {
8229: MatCheckPreallocated(mat,1);
8230: if (nullsp) {PetscObjectReference((PetscObject)nullsp);}
8231: MatNullSpaceDestroy(&mat->nearnullsp);
8232: mat->nearnullsp = nullsp;
8233: return(0);
8234: }
8236: /*@
8237: MatGetNearNullSpace -Get null space attached with MatSetNearNullSpace()
8239: Not Collective
8241: Input Parameters:
8242: . mat - the matrix
8244: Output Parameters:
8245: . nullsp - the null space object, NULL if not set
8247: Level: developer
8249: Concepts: null space^attaching to matrix
8251: .seealso: MatSetNearNullSpace(), MatGetNullSpace(), MatNullSpaceCreate()
8252: @*/
8253: PetscErrorCode MatGetNearNullSpace(Mat mat,MatNullSpace *nullsp)
8254: {
8259: MatCheckPreallocated(mat,1);
8260: *nullsp = mat->nearnullsp;
8261: return(0);
8262: }
8264: /*@C
8265: MatICCFactor - Performs in-place incomplete Cholesky factorization of matrix.
8267: Collective on Mat
8269: Input Parameters:
8270: + mat - the matrix
8271: . row - row/column permutation
8272: . fill - expected fill factor >= 1.0
8273: - level - level of fill, for ICC(k)
8275: Notes:
8276: Probably really in-place only when level of fill is zero, otherwise allocates
8277: new space to store factored matrix and deletes previous memory.
8279: Most users should employ the simplified KSP interface for linear solvers
8280: instead of working directly with matrix algebra routines such as this.
8281: See, e.g., KSPCreate().
8283: Level: developer
8285: Concepts: matrices^incomplete Cholesky factorization
8286: Concepts: Cholesky factorization
8288: .seealso: MatICCFactorSymbolic(), MatLUFactorNumeric(), MatCholeskyFactor()
8290: Developer Note: fortran interface is not autogenerated as the f90
8291: interface defintion cannot be generated correctly [due to MatFactorInfo]
8293: @*/
8294: PetscErrorCode MatICCFactor(Mat mat,IS row,const MatFactorInfo *info)
8295: {
8303: if (mat->rmap->N != mat->cmap->N) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONG,"matrix must be square");
8304: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
8305: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
8306: if (!mat->ops->iccfactor) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
8307: MatCheckPreallocated(mat,1);
8308: (*mat->ops->iccfactor)(mat,row,info);
8309: PetscObjectStateIncrease((PetscObject)mat);
8310: return(0);
8311: }
8313: /*@
8314: MatDiagonalScaleLocal - Scales columns of a matrix given the scaling values including the
8315: ghosted ones.
8317: Not Collective
8319: Input Parameters:
8320: + mat - the matrix
8321: - diag = the diagonal values, including ghost ones
8323: Level: developer
8325: Notes: Works only for MPIAIJ and MPIBAIJ matrices
8327: .seealso: MatDiagonalScale()
8328: @*/
8329: PetscErrorCode MatDiagonalScaleLocal(Mat mat,Vec diag)
8330: {
8332: PetscMPIInt size;
8339: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Matrix must be already assembled");
8340: PetscLogEventBegin(MAT_Scale,mat,0,0,0);
8341: MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);
8342: if (size == 1) {
8343: PetscInt n,m;
8344: VecGetSize(diag,&n);
8345: MatGetSize(mat,0,&m);
8346: if (m == n) {
8347: MatDiagonalScale(mat,0,diag);
8348: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only supported for sequential matrices when no ghost points/periodic conditions");
8349: } else {
8350: PetscUseMethod(mat,"MatDiagonalScaleLocal_C",(Mat,Vec),(mat,diag));
8351: }
8352: PetscLogEventEnd(MAT_Scale,mat,0,0,0);
8353: PetscObjectStateIncrease((PetscObject)mat);
8354: return(0);
8355: }
8357: /*@
8358: MatGetInertia - Gets the inertia from a factored matrix
8360: Collective on Mat
8362: Input Parameter:
8363: . mat - the matrix
8365: Output Parameters:
8366: + nneg - number of negative eigenvalues
8367: . nzero - number of zero eigenvalues
8368: - npos - number of positive eigenvalues
8370: Level: advanced
8372: Notes: Matrix must have been factored by MatCholeskyFactor()
8375: @*/
8376: PetscErrorCode MatGetInertia(Mat mat,PetscInt *nneg,PetscInt *nzero,PetscInt *npos)
8377: {
8383: if (!mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Unfactored matrix");
8384: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Numeric factor mat is not assembled");
8385: if (!mat->ops->getinertia) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
8386: (*mat->ops->getinertia)(mat,nneg,nzero,npos);
8387: return(0);
8388: }
8390: /* ----------------------------------------------------------------*/
8391: /*@C
8392: MatSolves - Solves A x = b, given a factored matrix, for a collection of vectors
8394: Neighbor-wise Collective on Mat and Vecs
8396: Input Parameters:
8397: + mat - the factored matrix
8398: - b - the right-hand-side vectors
8400: Output Parameter:
8401: . x - the result vectors
8403: Notes:
8404: The vectors b and x cannot be the same. I.e., one cannot
8405: call MatSolves(A,x,x).
8407: Notes:
8408: Most users should employ the simplified KSP interface for linear solvers
8409: instead of working directly with matrix algebra routines such as this.
8410: See, e.g., KSPCreate().
8412: Level: developer
8414: Concepts: matrices^triangular solves
8416: .seealso: MatSolveAdd(), MatSolveTranspose(), MatSolveTransposeAdd(), MatSolve()
8417: @*/
8418: PetscErrorCode MatSolves(Mat mat,Vecs b,Vecs x)
8419: {
8425: if (x == b) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_IDN,"x and b must be different vectors");
8426: if (!mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Unfactored matrix");
8427: if (!mat->rmap->N && !mat->cmap->N) return(0);
8429: if (!mat->ops->solves) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
8430: MatCheckPreallocated(mat,1);
8431: PetscLogEventBegin(MAT_Solves,mat,0,0,0);
8432: (*mat->ops->solves)(mat,b,x);
8433: PetscLogEventEnd(MAT_Solves,mat,0,0,0);
8434: return(0);
8435: }
8437: /*@
8438: MatIsSymmetric - Test whether a matrix is symmetric
8440: Collective on Mat
8442: Input Parameter:
8443: + A - the matrix to test
8444: - tol - difference between value and its transpose less than this amount counts as equal (use 0.0 for exact transpose)
8446: Output Parameters:
8447: . flg - the result
8449: Notes: For real numbers MatIsSymmetric() and MatIsHermitian() return identical results
8451: Level: intermediate
8453: Concepts: matrix^symmetry
8455: .seealso: MatTranspose(), MatIsTranspose(), MatIsHermitian(), MatIsStructurallySymmetric(), MatSetOption(), MatIsSymmetricKnown()
8456: @*/
8457: PetscErrorCode MatIsSymmetric(Mat A,PetscReal tol,PetscBool *flg)
8458: {
8465: if (!A->symmetric_set) {
8466: if (!A->ops->issymmetric) {
8467: MatType mattype;
8468: MatGetType(A,&mattype);
8469: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix of type <%s> does not support checking for symmetric",mattype);
8470: }
8471: (*A->ops->issymmetric)(A,tol,flg);
8472: if (!tol) {
8473: A->symmetric_set = PETSC_TRUE;
8474: A->symmetric = *flg;
8475: if (A->symmetric) {
8476: A->structurally_symmetric_set = PETSC_TRUE;
8477: A->structurally_symmetric = PETSC_TRUE;
8478: }
8479: }
8480: } else if (A->symmetric) {
8481: *flg = PETSC_TRUE;
8482: } else if (!tol) {
8483: *flg = PETSC_FALSE;
8484: } else {
8485: if (!A->ops->issymmetric) {
8486: MatType mattype;
8487: MatGetType(A,&mattype);
8488: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix of type <%s> does not support checking for symmetric",mattype);
8489: }
8490: (*A->ops->issymmetric)(A,tol,flg);
8491: }
8492: return(0);
8493: }
8495: /*@
8496: MatIsHermitian - Test whether a matrix is Hermitian
8498: Collective on Mat
8500: Input Parameter:
8501: + A - the matrix to test
8502: - tol - difference between value and its transpose less than this amount counts as equal (use 0.0 for exact Hermitian)
8504: Output Parameters:
8505: . flg - the result
8507: Level: intermediate
8509: Concepts: matrix^symmetry
8511: .seealso: MatTranspose(), MatIsTranspose(), MatIsHermitian(), MatIsStructurallySymmetric(), MatSetOption(),
8512: MatIsSymmetricKnown(), MatIsSymmetric()
8513: @*/
8514: PetscErrorCode MatIsHermitian(Mat A,PetscReal tol,PetscBool *flg)
8515: {
8522: if (!A->hermitian_set) {
8523: if (!A->ops->ishermitian) {
8524: MatType mattype;
8525: MatGetType(A,&mattype);
8526: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix of type <%s> does not support checking for hermitian",mattype);
8527: }
8528: (*A->ops->ishermitian)(A,tol,flg);
8529: if (!tol) {
8530: A->hermitian_set = PETSC_TRUE;
8531: A->hermitian = *flg;
8532: if (A->hermitian) {
8533: A->structurally_symmetric_set = PETSC_TRUE;
8534: A->structurally_symmetric = PETSC_TRUE;
8535: }
8536: }
8537: } else if (A->hermitian) {
8538: *flg = PETSC_TRUE;
8539: } else if (!tol) {
8540: *flg = PETSC_FALSE;
8541: } else {
8542: if (!A->ops->ishermitian) {
8543: MatType mattype;
8544: MatGetType(A,&mattype);
8545: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix of type <%s> does not support checking for hermitian",mattype);
8546: }
8547: (*A->ops->ishermitian)(A,tol,flg);
8548: }
8549: return(0);
8550: }
8552: /*@
8553: MatIsSymmetricKnown - Checks the flag on the matrix to see if it is symmetric.
8555: Not Collective
8557: Input Parameter:
8558: . A - the matrix to check
8560: Output Parameters:
8561: + set - if the symmetric flag is set (this tells you if the next flag is valid)
8562: - flg - the result
8564: Level: advanced
8566: Concepts: matrix^symmetry
8568: Note: Does not check the matrix values directly, so this may return unknown (set = PETSC_FALSE). Use MatIsSymmetric()
8569: if you want it explicitly checked
8571: .seealso: MatTranspose(), MatIsTranspose(), MatIsHermitian(), MatIsStructurallySymmetric(), MatSetOption(), MatIsSymmetric()
8572: @*/
8573: PetscErrorCode MatIsSymmetricKnown(Mat A,PetscBool *set,PetscBool *flg)
8574: {
8579: if (A->symmetric_set) {
8580: *set = PETSC_TRUE;
8581: *flg = A->symmetric;
8582: } else {
8583: *set = PETSC_FALSE;
8584: }
8585: return(0);
8586: }
8588: /*@
8589: MatIsHermitianKnown - Checks the flag on the matrix to see if it is hermitian.
8591: Not Collective
8593: Input Parameter:
8594: . A - the matrix to check
8596: Output Parameters:
8597: + set - if the hermitian flag is set (this tells you if the next flag is valid)
8598: - flg - the result
8600: Level: advanced
8602: Concepts: matrix^symmetry
8604: Note: Does not check the matrix values directly, so this may return unknown (set = PETSC_FALSE). Use MatIsHermitian()
8605: if you want it explicitly checked
8607: .seealso: MatTranspose(), MatIsTranspose(), MatIsHermitian(), MatIsStructurallySymmetric(), MatSetOption(), MatIsSymmetric()
8608: @*/
8609: PetscErrorCode MatIsHermitianKnown(Mat A,PetscBool *set,PetscBool *flg)
8610: {
8615: if (A->hermitian_set) {
8616: *set = PETSC_TRUE;
8617: *flg = A->hermitian;
8618: } else {
8619: *set = PETSC_FALSE;
8620: }
8621: return(0);
8622: }
8624: /*@
8625: MatIsStructurallySymmetric - Test whether a matrix is structurally symmetric
8627: Collective on Mat
8629: Input Parameter:
8630: . A - the matrix to test
8632: Output Parameters:
8633: . flg - the result
8635: Level: intermediate
8637: Concepts: matrix^symmetry
8639: .seealso: MatTranspose(), MatIsTranspose(), MatIsHermitian(), MatIsSymmetric(), MatSetOption()
8640: @*/
8641: PetscErrorCode MatIsStructurallySymmetric(Mat A,PetscBool *flg)
8642: {
8648: if (!A->structurally_symmetric_set) {
8649: if (!A->ops->isstructurallysymmetric) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Matrix does not support checking for structural symmetric");
8650: (*A->ops->isstructurallysymmetric)(A,&A->structurally_symmetric);
8652: A->structurally_symmetric_set = PETSC_TRUE;
8653: }
8654: *flg = A->structurally_symmetric;
8655: return(0);
8656: }
8658: /*@
8659: MatStashGetInfo - Gets how many values are currently in the matrix stash, i.e. need
8660: to be communicated to other processors during the MatAssemblyBegin/End() process
8662: Not collective
8664: Input Parameter:
8665: . vec - the vector
8667: Output Parameters:
8668: + nstash - the size of the stash
8669: . reallocs - the number of additional mallocs incurred.
8670: . bnstash - the size of the block stash
8671: - breallocs - the number of additional mallocs incurred.in the block stash
8673: Level: advanced
8675: .seealso: MatAssemblyBegin(), MatAssemblyEnd(), Mat, MatStashSetInitialSize()
8677: @*/
8678: PetscErrorCode MatStashGetInfo(Mat mat,PetscInt *nstash,PetscInt *reallocs,PetscInt *bnstash,PetscInt *breallocs)
8679: {
8683: MatStashGetInfo_Private(&mat->stash,nstash,reallocs);
8684: MatStashGetInfo_Private(&mat->bstash,bnstash,breallocs);
8685: return(0);
8686: }
8688: /*@C
8689: MatCreateVecs - Get vector(s) compatible with the matrix, i.e. with the same
8690: parallel layout
8692: Collective on Mat
8694: Input Parameter:
8695: . mat - the matrix
8697: Output Parameter:
8698: + right - (optional) vector that the matrix can be multiplied against
8699: - left - (optional) vector that the matrix vector product can be stored in
8701: Notes:
8702: The blocksize of the returned vectors is determined by the row and column block sizes set with MatSetBlockSizes() or the single blocksize (same for both) set by MatSetBlockSize().
8704: Notes: These are new vectors which are not owned by the Mat, they should be destroyed in VecDestroy() when no longer needed
8706: Level: advanced
8708: .seealso: MatCreate(), VecDestroy()
8709: @*/
8710: PetscErrorCode MatCreateVecs(Mat mat,Vec *right,Vec *left)
8711: {
8717: if (mat->ops->getvecs) {
8718: (*mat->ops->getvecs)(mat,right,left);
8719: } else {
8720: PetscInt rbs,cbs;
8721: MatGetBlockSizes(mat,&rbs,&cbs);
8722: if (right) {
8723: if (mat->cmap->n < 0) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"PetscLayout for columns not yet setup");
8724: VecCreate(PetscObjectComm((PetscObject)mat),right);
8725: VecSetSizes(*right,mat->cmap->n,PETSC_DETERMINE);
8726: VecSetBlockSize(*right,cbs);
8727: VecSetType(*right,VECSTANDARD);
8728: PetscLayoutReference(mat->cmap,&(*right)->map);
8729: }
8730: if (left) {
8731: if (mat->rmap->n < 0) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"PetscLayout for rows not yet setup");
8732: VecCreate(PetscObjectComm((PetscObject)mat),left);
8733: VecSetSizes(*left,mat->rmap->n,PETSC_DETERMINE);
8734: VecSetBlockSize(*left,rbs);
8735: VecSetType(*left,VECSTANDARD);
8736: PetscLayoutReference(mat->rmap,&(*left)->map);
8737: }
8738: }
8739: return(0);
8740: }
8742: /*@C
8743: MatFactorInfoInitialize - Initializes a MatFactorInfo data structure
8744: with default values.
8746: Not Collective
8748: Input Parameters:
8749: . info - the MatFactorInfo data structure
8752: Notes: The solvers are generally used through the KSP and PC objects, for example
8753: PCLU, PCILU, PCCHOLESKY, PCICC
8755: Level: developer
8757: .seealso: MatFactorInfo
8759: Developer Note: fortran interface is not autogenerated as the f90
8760: interface defintion cannot be generated correctly [due to MatFactorInfo]
8762: @*/
8764: PetscErrorCode MatFactorInfoInitialize(MatFactorInfo *info)
8765: {
8769: PetscMemzero(info,sizeof(MatFactorInfo));
8770: return(0);
8771: }
8773: /*@
8774: MatFactorSetSchurIS - Set indices corresponding to the Schur complement you wish to have computed
8776: Collective on Mat
8778: Input Parameters:
8779: + mat - the factored matrix
8780: - is - the index set defining the Schur indices (0-based)
8782: Notes: Call MatFactorSolveSchurComplement() or MatFactorSolveSchurComplementTranspose() after this call to solve a Schur complement system.
8784: You can call MatFactorGetSchurComplement() or MatFactorCreateSchurComplement() after this call.
8786: Level: developer
8788: Concepts:
8790: .seealso: MatGetFactor(), MatFactorGetSchurComplement(), MatFactorRestoreSchurComplement(), MatFactorCreateSchurComplement(), MatFactorSolveSchurComplement(),
8791: MatFactorSolveSchurComplementTranspose(), MatFactorSolveSchurComplement()
8793: @*/
8794: PetscErrorCode MatFactorSetSchurIS(Mat mat,IS is)
8795: {
8796: PetscErrorCode ierr,(*f)(Mat,IS);
8804: if (!mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
8805: PetscObjectQueryFunction((PetscObject)mat,"MatFactorSetSchurIS_C",&f);
8806: if (!f) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"The selected MatSolverType does not support Schur complement computation. You should use MATSOLVERMUMPS or MATSOLVERMKL_PARDISO");
8807: if (mat->schur) {
8808: MatDestroy(&mat->schur);
8809: }
8810: (*f)(mat,is);
8811: if (!mat->schur) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_PLIB,"Schur complement has not been created");
8812: MatFactorSetUpInPlaceSchur_Private(mat);
8813: return(0);
8814: }
8816: /*@
8817: MatFactorCreateSchurComplement - Create a Schur complement matrix object using Schur data computed during the factorization step
8819: Logically Collective on Mat
8821: Input Parameters:
8822: + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
8823: . S - location where to return the Schur complement, can be NULL
8824: - status - the status of the Schur complement matrix, can be NULL
8826: Notes:
8827: You must call MatFactorSetSchurIS() before calling this routine.
8829: The routine provides a copy of the Schur matrix stored within the solver data structures.
8830: The caller must destroy the object when it is no longer needed.
8831: If MatFactorInvertSchurComplement() has been called, the routine gets back the inverse.
8833: Use MatFactorGetSchurComplement() to get access to the Schur complement matrix inside the factored matrix instead of making a copy of it (which this function does)
8835: Developer Notes: The reason this routine exists is because the representation of the Schur complement within the factor matrix may be different than a standard PETSc
8836: matrix representation and we normally do not want to use the time or memory to make a copy as a regular PETSc matrix.
8838: See MatCreateSchurComplement() or MatGetSchurComplement() for ways to create virtual or approximate Schur complements.
8840: Level: advanced
8842: References:
8844: .seealso: MatGetFactor(), MatFactorSetSchurIS(), MatFactorGetSchurComplement(), MatFactorSchurStatus
8845: @*/
8846: PetscErrorCode MatFactorCreateSchurComplement(Mat F,Mat* S,MatFactorSchurStatus* status)
8847: {
8854: if (S) {
8855: PetscErrorCode (*f)(Mat,Mat*);
8857: PetscObjectQueryFunction((PetscObject)F,"MatFactorCreateSchurComplement_C",&f);
8858: if (f) {
8859: (*f)(F,S);
8860: } else {
8861: MatDuplicate(F->schur,MAT_COPY_VALUES,S);
8862: }
8863: }
8864: if (status) *status = F->schur_status;
8865: return(0);
8866: }
8868: /*@
8869: MatFactorGetSchurComplement - Gets access to a Schur complement matrix using the current Schur data within a factored matrix
8871: Logically Collective on Mat
8873: Input Parameters:
8874: + F - the factored matrix obtained by calling MatGetFactor()
8875: . *S - location where to return the Schur complement, can be NULL
8876: - status - the status of the Schur complement matrix, can be NULL
8878: Notes:
8879: You must call MatFactorSetSchurIS() before calling this routine.
8881: Schur complement mode is currently implemented for sequential matrices.
8882: The routine returns a the Schur Complement stored within the data strutures of the solver.
8883: If MatFactorInvertSchurComplement() has previously been called, the returned matrix is actually the inverse of the Schur complement.
8884: The returned matrix should not be destroyed; the caller should call MatFactorRestoreSchurComplement() when the object is no longer needed.
8886: Use MatFactorCreateSchurComplement() to create a copy of the Schur complement matrix that is within a factored matrix
8888: See MatCreateSchurComplement() or MatGetSchurComplement() for ways to create virtual or approximate Schur complements.
8890: Level: advanced
8892: References:
8894: .seealso: MatGetFactor(), MatFactorSetSchurIS(), MatFactorRestoreSchurComplement(), MatFactorCreateSchurComplement(), MatFactorSchurStatus
8895: @*/
8896: PetscErrorCode MatFactorGetSchurComplement(Mat F,Mat* S,MatFactorSchurStatus* status)
8897: {
8902: if (S) *S = F->schur;
8903: if (status) *status = F->schur_status;
8904: return(0);
8905: }
8907: /*@
8908: MatFactorRestoreSchurComplement - Restore the Schur complement matrix object obtained from a call to MatFactorGetSchurComplement
8910: Logically Collective on Mat
8912: Input Parameters:
8913: + F - the factored matrix obtained by calling MatGetFactor()
8914: . *S - location where the Schur complement is stored
8915: - status - the status of the Schur complement matrix (see MatFactorSchurStatus)
8917: Notes:
8919: Level: advanced
8921: References:
8923: .seealso: MatGetFactor(), MatFactorSetSchurIS(), MatFactorRestoreSchurComplement(), MatFactorCreateSchurComplement(), MatFactorSchurStatus
8924: @*/
8925: PetscErrorCode MatFactorRestoreSchurComplement(Mat F,Mat* S,MatFactorSchurStatus status)
8926: {
8931: if (S) {
8933: *S = NULL;
8934: }
8935: F->schur_status = status;
8936: MatFactorUpdateSchurStatus_Private(F);
8937: return(0);
8938: }
8940: /*@
8941: MatFactorSolveSchurComplementTranspose - Solve the transpose of the Schur complement system computed during the factorization step
8943: Logically Collective on Mat
8945: Input Parameters:
8946: + F - the factored matrix obtained by calling MatGetFactor()
8947: . rhs - location where the right hand side of the Schur complement system is stored
8948: - sol - location where the solution of the Schur complement system has to be returned
8950: Notes:
8951: The sizes of the vectors should match the size of the Schur complement
8953: Must be called after MatFactorSetSchurIS()
8955: Level: advanced
8957: References:
8959: .seealso: MatGetFactor(), MatFactorSetSchurIS(), MatFactorSolveSchurComplement()
8960: @*/
8961: PetscErrorCode MatFactorSolveSchurComplementTranspose(Mat F, Vec rhs, Vec sol)
8962: {
8974: MatFactorFactorizeSchurComplement(F);
8975: switch (F->schur_status) {
8976: case MAT_FACTOR_SCHUR_FACTORED:
8977: MatSolveTranspose(F->schur,rhs,sol);
8978: break;
8979: case MAT_FACTOR_SCHUR_INVERTED:
8980: MatMultTranspose(F->schur,rhs,sol);
8981: break;
8982: default:
8983: SETERRQ1(PetscObjectComm((PetscObject)F),PETSC_ERR_SUP,"Unhandled MatFactorSchurStatus %D",F->schur_status);
8984: break;
8985: }
8986: return(0);
8987: }
8989: /*@
8990: MatFactorSolveSchurComplement - Solve the Schur complement system computed during the factorization step
8992: Logically Collective on Mat
8994: Input Parameters:
8995: + F - the factored matrix obtained by calling MatGetFactor()
8996: . rhs - location where the right hand side of the Schur complement system is stored
8997: - sol - location where the solution of the Schur complement system has to be returned
8999: Notes:
9000: The sizes of the vectors should match the size of the Schur complement
9002: Must be called after MatFactorSetSchurIS()
9004: Level: advanced
9006: References:
9008: .seealso: MatGetFactor(), MatFactorSetSchurIS(), MatFactorSolveSchurComplementTranspose()
9009: @*/
9010: PetscErrorCode MatFactorSolveSchurComplement(Mat F, Vec rhs, Vec sol)
9011: {
9023: MatFactorFactorizeSchurComplement(F);
9024: switch (F->schur_status) {
9025: case MAT_FACTOR_SCHUR_FACTORED:
9026: MatSolve(F->schur,rhs,sol);
9027: break;
9028: case MAT_FACTOR_SCHUR_INVERTED:
9029: MatMult(F->schur,rhs,sol);
9030: break;
9031: default:
9032: SETERRQ1(PetscObjectComm((PetscObject)F),PETSC_ERR_SUP,"Unhandled MatFactorSchurStatus %D",F->schur_status);
9033: break;
9034: }
9035: return(0);
9036: }
9038: /*@
9039: MatFactorInvertSchurComplement - Invert the Schur complement matrix computed during the factorization step
9041: Logically Collective on Mat
9043: Input Parameters:
9044: + F - the factored matrix obtained by calling MatGetFactor()
9046: Notes: Must be called after MatFactorSetSchurIS().
9048: Call MatFactorGetSchurComplement() or MatFactorCreateSchurComplement() AFTER this call to actually compute the inverse and get access to it.
9050: Level: advanced
9052: References:
9054: .seealso: MatGetFactor(), MatFactorSetSchurIS(), MatFactorGetSchurComplement(), MatFactorCreateSchurComplement()
9055: @*/
9056: PetscErrorCode MatFactorInvertSchurComplement(Mat F)
9057: {
9063: if (F->schur_status == MAT_FACTOR_SCHUR_INVERTED) return(0);
9064: MatFactorFactorizeSchurComplement(F);
9065: MatFactorInvertSchurComplement_Private(F);
9066: F->schur_status = MAT_FACTOR_SCHUR_INVERTED;
9067: return(0);
9068: }
9070: /*@
9071: MatFactorFactorizeSchurComplement - Factorize the Schur complement matrix computed during the factorization step
9073: Logically Collective on Mat
9075: Input Parameters:
9076: + F - the factored matrix obtained by calling MatGetFactor()
9078: Notes: Must be called after MatFactorSetSchurIS().
9080: Level: advanced
9082: References:
9084: .seealso: MatGetFactor(), MatFactorSetSchurIS(), MatFactorInvertSchurComplement()
9085: @*/
9086: PetscErrorCode MatFactorFactorizeSchurComplement(Mat F)
9087: {
9093: if (F->schur_status == MAT_FACTOR_SCHUR_INVERTED || F->schur_status == MAT_FACTOR_SCHUR_FACTORED) return(0);
9094: MatFactorFactorizeSchurComplement_Private(F);
9095: F->schur_status = MAT_FACTOR_SCHUR_FACTORED;
9096: return(0);
9097: }
9099: /*@
9100: MatPtAP - Creates the matrix product C = P^T * A * P
9102: Neighbor-wise Collective on Mat
9104: Input Parameters:
9105: + A - the matrix
9106: . P - the projection matrix
9107: . scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
9108: - fill - expected fill as ratio of nnz(C)/(nnz(A) + nnz(P)), use PETSC_DEFAULT if you do not have a good estimate
9109: if the result is a dense matrix this is irrelevent
9111: Output Parameters:
9112: . C - the product matrix
9114: Notes:
9115: C will be created and must be destroyed by the user with MatDestroy().
9117: This routine is currently only implemented for pairs of sequential dense matrices, AIJ matrices and classes
9118: which inherit from AIJ.
9120: Level: intermediate
9122: .seealso: MatPtAPSymbolic(), MatPtAPNumeric(), MatMatMult(), MatRARt()
9123: @*/
9124: PetscErrorCode MatPtAP(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
9125: {
9127: PetscErrorCode (*fA)(Mat,Mat,MatReuse,PetscReal,Mat*);
9128: PetscErrorCode (*fP)(Mat,Mat,MatReuse,PetscReal,Mat*);
9129: PetscErrorCode (*ptap)(Mat,Mat,MatReuse,PetscReal,Mat*)=NULL;
9134: MatCheckPreallocated(A,1);
9135: if (scall == MAT_INPLACE_MATRIX) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Inplace product not supported");
9136: if (!A->assembled) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
9137: if (A->factortype) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
9140: MatCheckPreallocated(P,2);
9141: if (!P->assembled) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
9142: if (P->factortype) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
9144: if (A->rmap->N != A->cmap->N) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Matrix A must be square, %D != %D",A->rmap->N,A->cmap->N);
9145: if (P->rmap->N != A->cmap->N) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %D != %D",P->rmap->N,A->cmap->N);
9146: if (fill == PETSC_DEFAULT || fill == PETSC_DECIDE) fill = 2.0;
9147: if (fill < 1.0) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Expected fill=%g must be >= 1.0",(double)fill);
9149: if (scall == MAT_REUSE_MATRIX) {
9153: PetscLogEventBegin(MAT_PtAP,A,P,0,0);
9154: PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);
9155: (*(*C)->ops->ptapnumeric)(A,P,*C);
9156: PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);
9157: PetscLogEventEnd(MAT_PtAP,A,P,0,0);
9158: return(0);
9159: }
9161: if (fill == PETSC_DEFAULT || fill == PETSC_DECIDE) fill = 2.0;
9162: if (fill < 1.0) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Expected fill=%g must be >= 1.0",(double)fill);
9164: fA = A->ops->ptap;
9165: fP = P->ops->ptap;
9166: if (fP == fA) {
9167: if (!fA) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MatPtAP not supported for A of type %s",((PetscObject)A)->type_name);
9168: ptap = fA;
9169: } else {
9170: /* dispatch based on the type of A and P from their PetscObject's PetscFunctionLists. */
9171: char ptapname[256];
9172: PetscStrncpy(ptapname,"MatPtAP_",sizeof(ptapname));
9173: PetscStrlcat(ptapname,((PetscObject)A)->type_name,sizeof(ptapname));
9174: PetscStrlcat(ptapname,"_",sizeof(ptapname));
9175: PetscStrlcat(ptapname,((PetscObject)P)->type_name,sizeof(ptapname));
9176: PetscStrlcat(ptapname,"_C",sizeof(ptapname)); /* e.g., ptapname = "MatPtAP_seqdense_seqaij_C" */
9177: PetscObjectQueryFunction((PetscObject)P,ptapname,&ptap);
9178: if (!ptap) SETERRQ3(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"MatPtAP requires A, %s, to be compatible with P, %s (Misses composed function %s)",((PetscObject)A)->type_name,((PetscObject)P)->type_name,ptapname);
9179: }
9181: PetscLogEventBegin(MAT_PtAP,A,P,0,0);
9182: (*ptap)(A,P,scall,fill,C);
9183: PetscLogEventEnd(MAT_PtAP,A,P,0,0);
9184: return(0);
9185: }
9187: /*@
9188: MatPtAPNumeric - Computes the matrix product C = P^T * A * P
9190: Neighbor-wise Collective on Mat
9192: Input Parameters:
9193: + A - the matrix
9194: - P - the projection matrix
9196: Output Parameters:
9197: . C - the product matrix
9199: Notes:
9200: C must have been created by calling MatPtAPSymbolic and must be destroyed by
9201: the user using MatDeatroy().
9203: This routine is currently only implemented for pairs of AIJ matrices and classes
9204: which inherit from AIJ. C will be of type MATAIJ.
9206: Level: intermediate
9208: .seealso: MatPtAP(), MatPtAPSymbolic(), MatMatMultNumeric()
9209: @*/
9210: PetscErrorCode MatPtAPNumeric(Mat A,Mat P,Mat C)
9211: {
9217: if (!A->assembled) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
9218: if (A->factortype) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
9221: MatCheckPreallocated(P,2);
9222: if (!P->assembled) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
9223: if (P->factortype) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
9226: MatCheckPreallocated(C,3);
9227: if (C->factortype) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
9228: if (P->cmap->N!=C->rmap->N) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %D != %D",P->cmap->N,C->rmap->N);
9229: if (P->rmap->N!=A->cmap->N) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %D != %D",P->rmap->N,A->cmap->N);
9230: if (A->rmap->N!=A->cmap->N) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Matrix 'A' must be square, %D != %D",A->rmap->N,A->cmap->N);
9231: if (P->cmap->N!=C->cmap->N) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %D != %D",P->cmap->N,C->cmap->N);
9232: MatCheckPreallocated(A,1);
9234: if (!C->ops->ptapnumeric) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"MatPtAPNumeric implementation is missing. You should call MatPtAPSymbolic first");
9235: PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);
9236: (*C->ops->ptapnumeric)(A,P,C);
9237: PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);
9238: return(0);
9239: }
9241: /*@
9242: MatPtAPSymbolic - Creates the (i,j) structure of the matrix product C = P^T * A * P
9244: Neighbor-wise Collective on Mat
9246: Input Parameters:
9247: + A - the matrix
9248: - P - the projection matrix
9250: Output Parameters:
9251: . C - the (i,j) structure of the product matrix
9253: Notes:
9254: C will be created and must be destroyed by the user with MatDestroy().
9256: This routine is currently only implemented for pairs of SeqAIJ matrices and classes
9257: which inherit from SeqAIJ. C will be of type MATSEQAIJ. The product is computed using
9258: this (i,j) structure by calling MatPtAPNumeric().
9260: Level: intermediate
9262: .seealso: MatPtAP(), MatPtAPNumeric(), MatMatMultSymbolic()
9263: @*/
9264: PetscErrorCode MatPtAPSymbolic(Mat A,Mat P,PetscReal fill,Mat *C)
9265: {
9271: if (!A->assembled) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
9272: if (A->factortype) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
9273: if (fill <1.0) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Expected fill=%g must be >= 1.0",(double)fill);
9276: MatCheckPreallocated(P,2);
9277: if (!P->assembled) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
9278: if (P->factortype) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
9281: if (P->rmap->N!=A->cmap->N) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %D != %D",P->rmap->N,A->cmap->N);
9282: if (A->rmap->N!=A->cmap->N) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Matrix 'A' must be square, %D != %D",A->rmap->N,A->cmap->N);
9283: MatCheckPreallocated(A,1);
9285: if (!A->ops->ptapsymbolic) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MatType %s",((PetscObject)A)->type_name);
9286: PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);
9287: (*A->ops->ptapsymbolic)(A,P,fill,C);
9288: PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);
9290: /* MatSetBlockSize(*C,A->rmap->bs); NO! this is not always true -ma */
9291: return(0);
9292: }
9294: /*@
9295: MatRARt - Creates the matrix product C = R * A * R^T
9297: Neighbor-wise Collective on Mat
9299: Input Parameters:
9300: + A - the matrix
9301: . R - the projection matrix
9302: . scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
9303: - fill - expected fill as ratio of nnz(C)/nnz(A), use PETSC_DEFAULT if you do not have a good estimate
9304: if the result is a dense matrix this is irrelevent
9306: Output Parameters:
9307: . C - the product matrix
9309: Notes:
9310: C will be created and must be destroyed by the user with MatDestroy().
9312: This routine is currently only implemented for pairs of AIJ matrices and classes
9313: which inherit from AIJ. Due to PETSc sparse matrix block row distribution among processes,
9314: parallel MatRARt is implemented via explicit transpose of R, which could be very expensive.
9315: We recommend using MatPtAP().
9317: Level: intermediate
9319: .seealso: MatRARtSymbolic(), MatRARtNumeric(), MatMatMult(), MatPtAP()
9320: @*/
9321: PetscErrorCode MatRARt(Mat A,Mat R,MatReuse scall,PetscReal fill,Mat *C)
9322: {
9328: if (scall == MAT_INPLACE_MATRIX) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Inplace product not supported");
9329: if (!A->assembled) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
9330: if (A->factortype) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
9333: MatCheckPreallocated(R,2);
9334: if (!R->assembled) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
9335: if (R->factortype) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
9337: if (R->cmap->N!=A->rmap->N) SETERRQ2(PetscObjectComm((PetscObject)R),PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %D != %D",R->cmap->N,A->rmap->N);
9339: if (fill == PETSC_DEFAULT || fill == PETSC_DECIDE) fill = 2.0;
9340: if (fill < 1.0) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Expected fill=%g must be >= 1.0",(double)fill);
9341: MatCheckPreallocated(A,1);
9343: if (!A->ops->rart) {
9344: Mat Rt;
9345: MatTranspose(R,MAT_INITIAL_MATRIX,&Rt);
9346: MatMatMatMult(R,A,Rt,scall,fill,C);
9347: MatDestroy(&Rt);
9348: return(0);
9349: }
9350: PetscLogEventBegin(MAT_RARt,A,R,0,0);
9351: (*A->ops->rart)(A,R,scall,fill,C);
9352: PetscLogEventEnd(MAT_RARt,A,R,0,0);
9353: return(0);
9354: }
9356: /*@
9357: MatRARtNumeric - Computes the matrix product C = R * A * R^T
9359: Neighbor-wise Collective on Mat
9361: Input Parameters:
9362: + A - the matrix
9363: - R - the projection matrix
9365: Output Parameters:
9366: . C - the product matrix
9368: Notes:
9369: C must have been created by calling MatRARtSymbolic and must be destroyed by
9370: the user using MatDestroy().
9372: This routine is currently only implemented for pairs of AIJ matrices and classes
9373: which inherit from AIJ. C will be of type MATAIJ.
9375: Level: intermediate
9377: .seealso: MatRARt(), MatRARtSymbolic(), MatMatMultNumeric()
9378: @*/
9379: PetscErrorCode MatRARtNumeric(Mat A,Mat R,Mat C)
9380: {
9386: if (!A->assembled) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
9387: if (A->factortype) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
9390: MatCheckPreallocated(R,2);
9391: if (!R->assembled) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
9392: if (R->factortype) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
9395: MatCheckPreallocated(C,3);
9396: if (C->factortype) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
9397: if (R->rmap->N!=C->rmap->N) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %D != %D",R->rmap->N,C->rmap->N);
9398: if (R->cmap->N!=A->rmap->N) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %D != %D",R->cmap->N,A->rmap->N);
9399: if (A->rmap->N!=A->cmap->N) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Matrix 'A' must be square, %D != %D",A->rmap->N,A->cmap->N);
9400: if (R->rmap->N!=C->cmap->N) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %D != %D",R->rmap->N,C->cmap->N);
9401: MatCheckPreallocated(A,1);
9403: PetscLogEventBegin(MAT_RARtNumeric,A,R,0,0);
9404: (*A->ops->rartnumeric)(A,R,C);
9405: PetscLogEventEnd(MAT_RARtNumeric,A,R,0,0);
9406: return(0);
9407: }
9409: /*@
9410: MatRARtSymbolic - Creates the (i,j) structure of the matrix product C = R * A * R^T
9412: Neighbor-wise Collective on Mat
9414: Input Parameters:
9415: + A - the matrix
9416: - R - the projection matrix
9418: Output Parameters:
9419: . C - the (i,j) structure of the product matrix
9421: Notes:
9422: C will be created and must be destroyed by the user with MatDestroy().
9424: This routine is currently only implemented for pairs of SeqAIJ matrices and classes
9425: which inherit from SeqAIJ. C will be of type MATSEQAIJ. The product is computed using
9426: this (i,j) structure by calling MatRARtNumeric().
9428: Level: intermediate
9430: .seealso: MatRARt(), MatRARtNumeric(), MatMatMultSymbolic()
9431: @*/
9432: PetscErrorCode MatRARtSymbolic(Mat A,Mat R,PetscReal fill,Mat *C)
9433: {
9439: if (!A->assembled) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
9440: if (A->factortype) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
9441: if (fill <1.0) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Expected fill=%g must be >= 1.0",(double)fill);
9444: MatCheckPreallocated(R,2);
9445: if (!R->assembled) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
9446: if (R->factortype) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
9449: if (R->cmap->N!=A->rmap->N) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %D != %D",R->cmap->N,A->rmap->N);
9450: if (A->rmap->N!=A->cmap->N) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Matrix 'A' must be square, %D != %D",A->rmap->N,A->cmap->N);
9451: MatCheckPreallocated(A,1);
9452: PetscLogEventBegin(MAT_RARtSymbolic,A,R,0,0);
9453: (*A->ops->rartsymbolic)(A,R,fill,C);
9454: PetscLogEventEnd(MAT_RARtSymbolic,A,R,0,0);
9456: MatSetBlockSizes(*C,PetscAbs(R->rmap->bs),PetscAbs(R->rmap->bs));
9457: return(0);
9458: }
9460: /*@
9461: MatMatMult - Performs Matrix-Matrix Multiplication C=A*B.
9463: Neighbor-wise Collective on Mat
9465: Input Parameters:
9466: + A - the left matrix
9467: . B - the right matrix
9468: . scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
9469: - fill - expected fill as ratio of nnz(C)/(nnz(A) + nnz(B)), use PETSC_DEFAULT if you do not have a good estimate
9470: if the result is a dense matrix this is irrelevent
9472: Output Parameters:
9473: . C - the product matrix
9475: Notes:
9476: Unless scall is MAT_REUSE_MATRIX C will be created.
9478: MAT_REUSE_MATRIX can only be used if the matrices A and B have the same nonzero pattern as in the previous call and C was obtained from a previous
9479: call to this function with either MAT_INITIAL_MATRIX or MatMatMultSymbolic()
9481: To determine the correct fill value, run with -info and search for the string "Fill ratio" to see the value
9482: actually needed.
9484: If you have many matrices with the same non-zero structure to multiply, you
9485: should either
9486: $ 1) use MAT_REUSE_MATRIX in all calls but the first or
9487: $ 2) call MatMatMultSymbolic() once and then MatMatMultNumeric() for each product needed
9488: In the special case where matrix B (and hence C) are dense you can create the correctly sized matrix C yourself and then call this routine
9489: with MAT_REUSE_MATRIX, rather than first having MatMatMult() create it for you. You can NEVER do this if the matrix C is sparse.
9491: Level: intermediate
9493: .seealso: MatMatMultSymbolic(), MatMatMultNumeric(), MatTransposeMatMult(), MatMatTransposeMult(), MatPtAP()
9494: @*/
9495: PetscErrorCode MatMatMult(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
9496: {
9498: PetscErrorCode (*fA)(Mat,Mat,MatReuse,PetscReal,Mat*);
9499: PetscErrorCode (*fB)(Mat,Mat,MatReuse,PetscReal,Mat*);
9500: PetscErrorCode (*mult)(Mat,Mat,MatReuse,PetscReal,Mat*)=NULL;
9505: MatCheckPreallocated(A,1);
9506: if (!A->assembled) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
9507: if (A->factortype) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
9510: MatCheckPreallocated(B,2);
9511: if (!B->assembled) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
9512: if (B->factortype) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
9514: if (scall == MAT_INPLACE_MATRIX) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Inplace product not supported");
9515: if (B->rmap->N!=A->cmap->N) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %D != %D",B->rmap->N,A->cmap->N);
9516: if (scall == MAT_REUSE_MATRIX) {
9519: PetscLogEventBegin(MAT_MatMult,A,B,0,0);
9520: PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);
9521: (*(*C)->ops->matmultnumeric)(A,B,*C);
9522: PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);
9523: PetscLogEventEnd(MAT_MatMult,A,B,0,0);
9524: return(0);
9525: }
9526: if (fill == PETSC_DEFAULT || fill == PETSC_DECIDE) fill = 2.0;
9527: if (fill < 1.0) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Expected fill=%g must be >= 1.0",(double)fill);
9529: fA = A->ops->matmult;
9530: fB = B->ops->matmult;
9531: if (fB == fA) {
9532: if (!fB) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MatMatMult not supported for B of type %s",((PetscObject)B)->type_name);
9533: mult = fB;
9534: } else {
9535: /* dispatch based on the type of A and B from their PetscObject's PetscFunctionLists. */
9536: char multname[256];
9537: PetscStrncpy(multname,"MatMatMult_",sizeof(multname));
9538: PetscStrlcat(multname,((PetscObject)A)->type_name,sizeof(multname));
9539: PetscStrlcat(multname,"_",sizeof(multname));
9540: PetscStrlcat(multname,((PetscObject)B)->type_name,sizeof(multname));
9541: PetscStrlcat(multname,"_C",sizeof(multname)); /* e.g., multname = "MatMatMult_seqdense_seqaij_C" */
9542: PetscObjectQueryFunction((PetscObject)B,multname,&mult);
9543: if (!mult) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"MatMatMult requires A, %s, to be compatible with B, %s",((PetscObject)A)->type_name,((PetscObject)B)->type_name);
9544: }
9545: PetscLogEventBegin(MAT_MatMult,A,B,0,0);
9546: (*mult)(A,B,scall,fill,C);
9547: PetscLogEventEnd(MAT_MatMult,A,B,0,0);
9548: return(0);
9549: }
9551: /*@
9552: MatMatMultSymbolic - Performs construction, preallocation, and computes the ij structure
9553: of the matrix-matrix product C=A*B. Call this routine before calling MatMatMultNumeric().
9555: Neighbor-wise Collective on Mat
9557: Input Parameters:
9558: + A - the left matrix
9559: . B - the right matrix
9560: - fill - expected fill as ratio of nnz(C)/(nnz(A) + nnz(B)), use PETSC_DEFAULT if you do not have a good estimate,
9561: if C is a dense matrix this is irrelevent
9563: Output Parameters:
9564: . C - the product matrix
9566: Notes:
9567: Unless scall is MAT_REUSE_MATRIX C will be created.
9569: To determine the correct fill value, run with -info and search for the string "Fill ratio" to see the value
9570: actually needed.
9572: This routine is currently implemented for
9573: - pairs of AIJ matrices and classes which inherit from AIJ, C will be of type AIJ
9574: - pairs of AIJ (A) and Dense (B) matrix, C will be of type Dense.
9575: - pairs of Dense (A) and AIJ (B) matrix, C will be of type Dense.
9577: Level: intermediate
9579: Developers Note: There are ways to estimate the number of nonzeros in the resulting product, see for example, http://arxiv.org/abs/1006.4173
9580: We should incorporate them into PETSc.
9582: .seealso: MatMatMult(), MatMatMultNumeric()
9583: @*/
9584: PetscErrorCode MatMatMultSymbolic(Mat A,Mat B,PetscReal fill,Mat *C)
9585: {
9587: PetscErrorCode (*Asymbolic)(Mat,Mat,PetscReal,Mat*);
9588: PetscErrorCode (*Bsymbolic)(Mat,Mat,PetscReal,Mat*);
9589: PetscErrorCode (*symbolic)(Mat,Mat,PetscReal,Mat*)=NULL;
9594: if (!A->assembled) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
9595: if (A->factortype) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
9599: MatCheckPreallocated(B,2);
9600: if (!B->assembled) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
9601: if (B->factortype) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
9604: if (B->rmap->N!=A->cmap->N) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %D != %D",B->rmap->N,A->cmap->N);
9605: if (fill == PETSC_DEFAULT) fill = 2.0;
9606: if (fill < 1.0) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Expected fill=%g must be > 1.0",(double)fill);
9607: MatCheckPreallocated(A,1);
9609: Asymbolic = A->ops->matmultsymbolic;
9610: Bsymbolic = B->ops->matmultsymbolic;
9611: if (Asymbolic == Bsymbolic) {
9612: if (!Bsymbolic) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"C=A*B not implemented for B of type %s",((PetscObject)B)->type_name);
9613: symbolic = Bsymbolic;
9614: } else { /* dispatch based on the type of A and B */
9615: char symbolicname[256];
9616: PetscStrncpy(symbolicname,"MatMatMultSymbolic_",sizeof(symbolicname));
9617: PetscStrlcat(symbolicname,((PetscObject)A)->type_name,sizeof(symbolicname));
9618: PetscStrlcat(symbolicname,"_",sizeof(symbolicname));
9619: PetscStrlcat(symbolicname,((PetscObject)B)->type_name,sizeof(symbolicname));
9620: PetscStrlcat(symbolicname,"_C",sizeof(symbolicname));
9621: PetscObjectQueryFunction((PetscObject)B,symbolicname,&symbolic);
9622: if (!symbolic) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"MatMatMultSymbolic requires A, %s, to be compatible with B, %s",((PetscObject)A)->type_name,((PetscObject)B)->type_name);
9623: }
9624: PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);
9625: (*symbolic)(A,B,fill,C);
9626: PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);
9627: return(0);
9628: }
9630: /*@
9631: MatMatMultNumeric - Performs the numeric matrix-matrix product.
9632: Call this routine after first calling MatMatMultSymbolic().
9634: Neighbor-wise Collective on Mat
9636: Input Parameters:
9637: + A - the left matrix
9638: - B - the right matrix
9640: Output Parameters:
9641: . C - the product matrix, which was created by from MatMatMultSymbolic() or a call to MatMatMult().
9643: Notes:
9644: C must have been created with MatMatMultSymbolic().
9646: This routine is currently implemented for
9647: - pairs of AIJ matrices and classes which inherit from AIJ, C will be of type MATAIJ.
9648: - pairs of AIJ (A) and Dense (B) matrix, C will be of type Dense.
9649: - pairs of Dense (A) and AIJ (B) matrix, C will be of type Dense.
9651: Level: intermediate
9653: .seealso: MatMatMult(), MatMatMultSymbolic()
9654: @*/
9655: PetscErrorCode MatMatMultNumeric(Mat A,Mat B,Mat C)
9656: {
9660: MatMatMult(A,B,MAT_REUSE_MATRIX,0.0,&C);
9661: return(0);
9662: }
9664: /*@
9665: MatMatTransposeMult - Performs Matrix-Matrix Multiplication C=A*B^T.
9667: Neighbor-wise Collective on Mat
9669: Input Parameters:
9670: + A - the left matrix
9671: . B - the right matrix
9672: . scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
9673: - fill - expected fill as ratio of nnz(C)/(nnz(A) + nnz(B)), use PETSC_DEFAULT if not known
9675: Output Parameters:
9676: . C - the product matrix
9678: Notes:
9679: C will be created if MAT_INITIAL_MATRIX and must be destroyed by the user with MatDestroy().
9681: MAT_REUSE_MATRIX can only be used if the matrices A and B have the same nonzero pattern as in the previous call
9683: To determine the correct fill value, run with -info and search for the string "Fill ratio" to see the value
9684: actually needed.
9686: This routine is currently only implemented for pairs of SeqAIJ matrices and for the SeqDense class.
9688: Level: intermediate
9690: .seealso: MatMatTransposeMultSymbolic(), MatMatTransposeMultNumeric(), MatMatMult(), MatTransposeMatMult() MatPtAP()
9691: @*/
9692: PetscErrorCode MatMatTransposeMult(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
9693: {
9695: PetscErrorCode (*fA)(Mat,Mat,MatReuse,PetscReal,Mat*);
9696: PetscErrorCode (*fB)(Mat,Mat,MatReuse,PetscReal,Mat*);
9701: if (scall == MAT_INPLACE_MATRIX) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Inplace product not supported");
9702: if (!A->assembled) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
9703: if (A->factortype) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
9706: MatCheckPreallocated(B,2);
9707: if (!B->assembled) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
9708: if (B->factortype) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
9710: if (B->cmap->N!=A->cmap->N) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, AN %D != BN %D",A->cmap->N,B->cmap->N);
9711: if (fill == PETSC_DEFAULT || fill == PETSC_DECIDE) fill = 2.0;
9712: if (fill < 1.0) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Expected fill=%g must be > 1.0",(double)fill);
9713: MatCheckPreallocated(A,1);
9715: fA = A->ops->mattransposemult;
9716: if (!fA) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MatMatTransposeMult not supported for A of type %s",((PetscObject)A)->type_name);
9717: fB = B->ops->mattransposemult;
9718: if (!fB) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MatMatTransposeMult not supported for B of type %s",((PetscObject)B)->type_name);
9719: if (fB!=fA) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"MatMatTransposeMult requires A, %s, to be compatible with B, %s",((PetscObject)A)->type_name,((PetscObject)B)->type_name);
9721: PetscLogEventBegin(MAT_MatTransposeMult,A,B,0,0);
9722: if (scall == MAT_INITIAL_MATRIX) {
9723: PetscLogEventBegin(MAT_MatTransposeMultSymbolic,A,B,0,0);
9724: (*A->ops->mattransposemultsymbolic)(A,B,fill,C);
9725: PetscLogEventEnd(MAT_MatTransposeMultSymbolic,A,B,0,0);
9726: }
9727: PetscLogEventBegin(MAT_MatTransposeMultNumeric,A,B,0,0);
9728: (*A->ops->mattransposemultnumeric)(A,B,*C);
9729: PetscLogEventEnd(MAT_MatTransposeMultNumeric,A,B,0,0);
9730: PetscLogEventEnd(MAT_MatTransposeMult,A,B,0,0);
9731: return(0);
9732: }
9734: /*@
9735: MatTransposeMatMult - Performs Matrix-Matrix Multiplication C=A^T*B.
9737: Neighbor-wise Collective on Mat
9739: Input Parameters:
9740: + A - the left matrix
9741: . B - the right matrix
9742: . scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
9743: - fill - expected fill as ratio of nnz(C)/(nnz(A) + nnz(B)), use PETSC_DEFAULT if not known
9745: Output Parameters:
9746: . C - the product matrix
9748: Notes:
9749: C will be created if MAT_INITIAL_MATRIX and must be destroyed by the user with MatDestroy().
9751: MAT_REUSE_MATRIX can only be used if the matrices A and B have the same nonzero pattern as in the previous call
9753: To determine the correct fill value, run with -info and search for the string "Fill ratio" to see the value
9754: actually needed.
9756: This routine is currently implemented for pairs of AIJ matrices and pairs of SeqDense matrices and classes
9757: which inherit from SeqAIJ. C will be of same type as the input matrices.
9759: Level: intermediate
9761: .seealso: MatTransposeMatMultSymbolic(), MatTransposeMatMultNumeric(), MatMatMult(), MatMatTransposeMult(), MatPtAP()
9762: @*/
9763: PetscErrorCode MatTransposeMatMult(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
9764: {
9766: PetscErrorCode (*fA)(Mat,Mat,MatReuse,PetscReal,Mat*);
9767: PetscErrorCode (*fB)(Mat,Mat,MatReuse,PetscReal,Mat*);
9768: PetscErrorCode (*transposematmult)(Mat,Mat,MatReuse,PetscReal,Mat*) = NULL;
9773: if (scall == MAT_INPLACE_MATRIX) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Inplace product not supported");
9774: if (!A->assembled) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
9775: if (A->factortype) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
9778: MatCheckPreallocated(B,2);
9779: if (!B->assembled) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
9780: if (B->factortype) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
9782: if (B->rmap->N!=A->rmap->N) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %D != %D",B->rmap->N,A->rmap->N);
9783: if (fill == PETSC_DEFAULT || fill == PETSC_DECIDE) fill = 2.0;
9784: if (fill < 1.0) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Expected fill=%g must be > 1.0",(double)fill);
9785: MatCheckPreallocated(A,1);
9787: fA = A->ops->transposematmult;
9788: fB = B->ops->transposematmult;
9789: if (fB==fA) {
9790: if (!fA) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MatTransposeMatMult not supported for A of type %s",((PetscObject)A)->type_name);
9791: transposematmult = fA;
9792: } else {
9793: /* dispatch based on the type of A and B from their PetscObject's PetscFunctionLists. */
9794: char multname[256];
9795: PetscStrncpy(multname,"MatTransposeMatMult_",sizeof(multname));
9796: PetscStrlcat(multname,((PetscObject)A)->type_name,sizeof(multname));
9797: PetscStrlcat(multname,"_",sizeof(multname));
9798: PetscStrlcat(multname,((PetscObject)B)->type_name,sizeof(multname));
9799: PetscStrlcat(multname,"_C",sizeof(multname)); /* e.g., multname = "MatMatMult_seqdense_seqaij_C" */
9800: PetscObjectQueryFunction((PetscObject)B,multname,&transposematmult);
9801: if (!transposematmult) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"MatTransposeMatMult requires A, %s, to be compatible with B, %s",((PetscObject)A)->type_name,((PetscObject)B)->type_name);
9802: }
9803: PetscLogEventBegin(MAT_TransposeMatMult,A,B,0,0);
9804: (*transposematmult)(A,B,scall,fill,C);
9805: PetscLogEventEnd(MAT_TransposeMatMult,A,B,0,0);
9806: return(0);
9807: }
9809: /*@
9810: MatMatMatMult - Performs Matrix-Matrix-Matrix Multiplication D=A*B*C.
9812: Neighbor-wise Collective on Mat
9814: Input Parameters:
9815: + A - the left matrix
9816: . B - the middle matrix
9817: . C - the right matrix
9818: . scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
9819: - fill - expected fill as ratio of nnz(D)/(nnz(A) + nnz(B)+nnz(C)), use PETSC_DEFAULT if you do not have a good estimate
9820: if the result is a dense matrix this is irrelevent
9822: Output Parameters:
9823: . D - the product matrix
9825: Notes:
9826: Unless scall is MAT_REUSE_MATRIX D will be created.
9828: MAT_REUSE_MATRIX can only be used if the matrices A, B and C have the same nonzero pattern as in the previous call
9830: To determine the correct fill value, run with -info and search for the string "Fill ratio" to see the value
9831: actually needed.
9833: If you have many matrices with the same non-zero structure to multiply, you
9834: should use MAT_REUSE_MATRIX in all calls but the first or
9836: Level: intermediate
9838: .seealso: MatMatMult, MatPtAP()
9839: @*/
9840: PetscErrorCode MatMatMatMult(Mat A,Mat B,Mat C,MatReuse scall,PetscReal fill,Mat *D)
9841: {
9843: PetscErrorCode (*fA)(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
9844: PetscErrorCode (*fB)(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
9845: PetscErrorCode (*fC)(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
9846: PetscErrorCode (*mult)(Mat,Mat,Mat,MatReuse,PetscReal,Mat*)=NULL;
9851: MatCheckPreallocated(A,1);
9852: if (scall == MAT_INPLACE_MATRIX) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Inplace product not supported");
9853: if (!A->assembled) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
9854: if (A->factortype) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
9857: MatCheckPreallocated(B,2);
9858: if (!B->assembled) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
9859: if (B->factortype) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
9862: MatCheckPreallocated(C,3);
9863: if (!C->assembled) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
9864: if (C->factortype) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
9865: if (B->rmap->N!=A->cmap->N) SETERRQ2(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %D != %D",B->rmap->N,A->cmap->N);
9866: if (C->rmap->N!=B->cmap->N) SETERRQ2(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %D != %D",C->rmap->N,B->cmap->N);
9867: if (scall == MAT_REUSE_MATRIX) {
9870: PetscLogEventBegin(MAT_MatMatMult,A,B,0,0);
9871: (*(*D)->ops->matmatmult)(A,B,C,scall,fill,D);
9872: PetscLogEventEnd(MAT_MatMatMult,A,B,0,0);
9873: return(0);
9874: }
9875: if (fill == PETSC_DEFAULT || fill == PETSC_DECIDE) fill = 2.0;
9876: if (fill < 1.0) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Expected fill=%g must be >= 1.0",(double)fill);
9878: fA = A->ops->matmatmult;
9879: fB = B->ops->matmatmult;
9880: fC = C->ops->matmatmult;
9881: if (fA == fB && fA == fC) {
9882: if (!fA) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MatMatMatMult not supported for A of type %s",((PetscObject)A)->type_name);
9883: mult = fA;
9884: } else {
9885: /* dispatch based on the type of A, B and C from their PetscObject's PetscFunctionLists. */
9886: char multname[256];
9887: PetscStrncpy(multname,"MatMatMatMult_",sizeof(multname));
9888: PetscStrlcat(multname,((PetscObject)A)->type_name,sizeof(multname));
9889: PetscStrlcat(multname,"_",sizeof(multname));
9890: PetscStrlcat(multname,((PetscObject)B)->type_name,sizeof(multname));
9891: PetscStrlcat(multname,"_",sizeof(multname));
9892: PetscStrlcat(multname,((PetscObject)C)->type_name,sizeof(multname));
9893: PetscStrlcat(multname,"_C",sizeof(multname));
9894: PetscObjectQueryFunction((PetscObject)B,multname,&mult);
9895: if (!mult) SETERRQ3(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"MatMatMatMult requires A, %s, to be compatible with B, %s, C, %s",((PetscObject)A)->type_name,((PetscObject)B)->type_name,((PetscObject)C)->type_name);
9896: }
9897: PetscLogEventBegin(MAT_MatMatMult,A,B,0,0);
9898: (*mult)(A,B,C,scall,fill,D);
9899: PetscLogEventEnd(MAT_MatMatMult,A,B,0,0);
9900: return(0);
9901: }
9903: /*@
9904: MatCreateRedundantMatrix - Create redundant matrices and put them into processors of subcommunicators.
9906: Collective on Mat
9908: Input Parameters:
9909: + mat - the matrix
9910: . nsubcomm - the number of subcommunicators (= number of redundant parallel or sequential matrices)
9911: . subcomm - MPI communicator split from the communicator where mat resides in (or MPI_COMM_NULL if nsubcomm is used)
9912: - reuse - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
9914: Output Parameter:
9915: . matredundant - redundant matrix
9917: Notes:
9918: MAT_REUSE_MATRIX can only be used when the nonzero structure of the
9919: original matrix has not changed from that last call to MatCreateRedundantMatrix().
9921: This routine creates the duplicated matrices in subcommunicators; you should NOT create them before
9922: calling it.
9924: Level: advanced
9926: Concepts: subcommunicator
9927: Concepts: duplicate matrix
9929: .seealso: MatDestroy()
9930: @*/
9931: PetscErrorCode MatCreateRedundantMatrix(Mat mat,PetscInt nsubcomm,MPI_Comm subcomm,MatReuse reuse,Mat *matredundant)
9932: {
9934: MPI_Comm comm;
9935: PetscMPIInt size;
9936: PetscInt mloc_sub,nloc_sub,rstart,rend,M=mat->rmap->N,N=mat->cmap->N,bs=mat->rmap->bs;
9937: Mat_Redundant *redund=NULL;
9938: PetscSubcomm psubcomm=NULL;
9939: MPI_Comm subcomm_in=subcomm;
9940: Mat *matseq;
9941: IS isrow,iscol;
9942: PetscBool newsubcomm=PETSC_FALSE;
9946: if (nsubcomm && reuse == MAT_REUSE_MATRIX) {
9949: }
9951: MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);
9952: if (size == 1 || nsubcomm == 1) {
9953: if (reuse == MAT_INITIAL_MATRIX) {
9954: MatDuplicate(mat,MAT_COPY_VALUES,matredundant);
9955: } else {
9956: if (*matredundant == mat) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"MAT_REUSE_MATRIX means reuse the matrix passed in as the final argument, not the original matrix");
9957: MatCopy(mat,*matredundant,SAME_NONZERO_PATTERN);
9958: }
9959: return(0);
9960: }
9962: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
9963: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
9964: MatCheckPreallocated(mat,1);
9966: PetscLogEventBegin(MAT_RedundantMat,mat,0,0,0);
9967: if (subcomm_in == MPI_COMM_NULL && reuse == MAT_INITIAL_MATRIX) { /* get subcomm if user does not provide subcomm */
9968: /* create psubcomm, then get subcomm */
9969: PetscObjectGetComm((PetscObject)mat,&comm);
9970: MPI_Comm_size(comm,&size);
9971: if (nsubcomm < 1 || nsubcomm > size) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"nsubcomm must between 1 and %D",size);
9973: PetscSubcommCreate(comm,&psubcomm);
9974: PetscSubcommSetNumber(psubcomm,nsubcomm);
9975: PetscSubcommSetType(psubcomm,PETSC_SUBCOMM_CONTIGUOUS);
9976: PetscSubcommSetFromOptions(psubcomm);
9977: PetscCommDuplicate(PetscSubcommChild(psubcomm),&subcomm,NULL);
9978: newsubcomm = PETSC_TRUE;
9979: PetscSubcommDestroy(&psubcomm);
9980: }
9982: /* get isrow, iscol and a local sequential matrix matseq[0] */
9983: if (reuse == MAT_INITIAL_MATRIX) {
9984: mloc_sub = PETSC_DECIDE;
9985: nloc_sub = PETSC_DECIDE;
9986: if (bs < 1) {
9987: PetscSplitOwnership(subcomm,&mloc_sub,&M);
9988: PetscSplitOwnership(subcomm,&nloc_sub,&N);
9989: } else {
9990: PetscSplitOwnershipBlock(subcomm,bs,&mloc_sub,&M);
9991: PetscSplitOwnershipBlock(subcomm,bs,&nloc_sub,&N);
9992: }
9993: MPI_Scan(&mloc_sub,&rend,1,MPIU_INT,MPI_SUM,subcomm);
9994: rstart = rend - mloc_sub;
9995: ISCreateStride(PETSC_COMM_SELF,mloc_sub,rstart,1,&isrow);
9996: ISCreateStride(PETSC_COMM_SELF,N,0,1,&iscol);
9997: } else { /* reuse == MAT_REUSE_MATRIX */
9998: if (*matredundant == mat) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"MAT_REUSE_MATRIX means reuse the matrix passed in as the final argument, not the original matrix");
9999: /* retrieve subcomm */
10000: PetscObjectGetComm((PetscObject)(*matredundant),&subcomm);
10001: redund = (*matredundant)->redundant;
10002: isrow = redund->isrow;
10003: iscol = redund->iscol;
10004: matseq = redund->matseq;
10005: }
10006: MatCreateSubMatrices(mat,1,&isrow,&iscol,reuse,&matseq);
10008: /* get matredundant over subcomm */
10009: if (reuse == MAT_INITIAL_MATRIX) {
10010: MatCreateMPIMatConcatenateSeqMat(subcomm,matseq[0],nloc_sub,reuse,matredundant);
10012: /* create a supporting struct and attach it to C for reuse */
10013: PetscNewLog(*matredundant,&redund);
10014: (*matredundant)->redundant = redund;
10015: redund->isrow = isrow;
10016: redund->iscol = iscol;
10017: redund->matseq = matseq;
10018: if (newsubcomm) {
10019: redund->subcomm = subcomm;
10020: } else {
10021: redund->subcomm = MPI_COMM_NULL;
10022: }
10023: } else {
10024: MatCreateMPIMatConcatenateSeqMat(subcomm,matseq[0],PETSC_DECIDE,reuse,matredundant);
10025: }
10026: PetscLogEventEnd(MAT_RedundantMat,mat,0,0,0);
10027: return(0);
10028: }
10030: /*@C
10031: MatGetMultiProcBlock - Create multiple [bjacobi] 'parallel submatrices' from
10032: a given 'mat' object. Each submatrix can span multiple procs.
10034: Collective on Mat
10036: Input Parameters:
10037: + mat - the matrix
10038: . subcomm - the subcommunicator obtained by com_split(comm)
10039: - scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
10041: Output Parameter:
10042: . subMat - 'parallel submatrices each spans a given subcomm
10044: Notes:
10045: The submatrix partition across processors is dictated by 'subComm' a
10046: communicator obtained by com_split(comm). The comm_split
10047: is not restriced to be grouped with consecutive original ranks.
10049: Due the comm_split() usage, the parallel layout of the submatrices
10050: map directly to the layout of the original matrix [wrt the local
10051: row,col partitioning]. So the original 'DiagonalMat' naturally maps
10052: into the 'DiagonalMat' of the subMat, hence it is used directly from
10053: the subMat. However the offDiagMat looses some columns - and this is
10054: reconstructed with MatSetValues()
10056: Level: advanced
10058: Concepts: subcommunicator
10059: Concepts: submatrices
10061: .seealso: MatCreateSubMatrices()
10062: @*/
10063: PetscErrorCode MatGetMultiProcBlock(Mat mat, MPI_Comm subComm, MatReuse scall,Mat *subMat)
10064: {
10066: PetscMPIInt commsize,subCommSize;
10069: MPI_Comm_size(PetscObjectComm((PetscObject)mat),&commsize);
10070: MPI_Comm_size(subComm,&subCommSize);
10071: if (subCommSize > commsize) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_OUTOFRANGE,"CommSize %D < SubCommZize %D",commsize,subCommSize);
10073: if (scall == MAT_REUSE_MATRIX && *subMat == mat) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"MAT_REUSE_MATRIX means reuse the matrix passed in as the final argument, not the original matrix");
10074: PetscLogEventBegin(MAT_GetMultiProcBlock,mat,0,0,0);
10075: (*mat->ops->getmultiprocblock)(mat,subComm,scall,subMat);
10076: PetscLogEventEnd(MAT_GetMultiProcBlock,mat,0,0,0);
10077: return(0);
10078: }
10080: /*@
10081: MatGetLocalSubMatrix - Gets a reference to a submatrix specified in local numbering
10083: Not Collective
10085: Input Arguments:
10086: mat - matrix to extract local submatrix from
10087: isrow - local row indices for submatrix
10088: iscol - local column indices for submatrix
10090: Output Arguments:
10091: submat - the submatrix
10093: Level: intermediate
10095: Notes:
10096: The submat should be returned with MatRestoreLocalSubMatrix().
10098: Depending on the format of mat, the returned submat may not implement MatMult(). Its communicator may be
10099: the same as mat, it may be PETSC_COMM_SELF, or some other subcomm of mat's.
10101: The submat always implements MatSetValuesLocal(). If isrow and iscol have the same block size, then
10102: MatSetValuesBlockedLocal() will also be implemented.
10104: The mat must have had a ISLocalToGlobalMapping provided to it with MatSetLocalToGlobalMapping(). Note that
10105: matrices obtained with DMCreateMat() generally already have the local to global mapping provided.
10107: .seealso: MatRestoreLocalSubMatrix(), MatCreateLocalRef(), MatSetLocalToGlobalMapping()
10108: @*/
10109: PetscErrorCode MatGetLocalSubMatrix(Mat mat,IS isrow,IS iscol,Mat *submat)
10110: {
10119: if (!mat->rmap->mapping) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Matrix must have local to global mapping provided before this call");
10121: if (mat->ops->getlocalsubmatrix) {
10122: (*mat->ops->getlocalsubmatrix)(mat,isrow,iscol,submat);
10123: } else {
10124: MatCreateLocalRef(mat,isrow,iscol,submat);
10125: }
10126: return(0);
10127: }
10129: /*@
10130: MatRestoreLocalSubMatrix - Restores a reference to a submatrix specified in local numbering
10132: Not Collective
10134: Input Arguments:
10135: mat - matrix to extract local submatrix from
10136: isrow - local row indices for submatrix
10137: iscol - local column indices for submatrix
10138: submat - the submatrix
10140: Level: intermediate
10142: .seealso: MatGetLocalSubMatrix()
10143: @*/
10144: PetscErrorCode MatRestoreLocalSubMatrix(Mat mat,IS isrow,IS iscol,Mat *submat)
10145: {
10154: if (*submat) {
10156: }
10158: if (mat->ops->restorelocalsubmatrix) {
10159: (*mat->ops->restorelocalsubmatrix)(mat,isrow,iscol,submat);
10160: } else {
10161: MatDestroy(submat);
10162: }
10163: *submat = NULL;
10164: return(0);
10165: }
10167: /* --------------------------------------------------------*/
10168: /*@
10169: MatFindZeroDiagonals - Finds all the rows of a matrix that have zero or no diagonal entry in the matrix
10171: Collective on Mat
10173: Input Parameter:
10174: . mat - the matrix
10176: Output Parameter:
10177: . is - if any rows have zero diagonals this contains the list of them
10179: Level: developer
10181: Concepts: matrix-vector product
10183: .seealso: MatMultTranspose(), MatMultAdd(), MatMultTransposeAdd()
10184: @*/
10185: PetscErrorCode MatFindZeroDiagonals(Mat mat,IS *is)
10186: {
10192: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
10193: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
10195: if (!mat->ops->findzerodiagonals) {
10196: Vec diag;
10197: const PetscScalar *a;
10198: PetscInt *rows;
10199: PetscInt rStart, rEnd, r, nrow = 0;
10201: MatCreateVecs(mat, &diag, NULL);
10202: MatGetDiagonal(mat, diag);
10203: MatGetOwnershipRange(mat, &rStart, &rEnd);
10204: VecGetArrayRead(diag, &a);
10205: for (r = 0; r < rEnd-rStart; ++r) if (a[r] == 0.0) ++nrow;
10206: PetscMalloc1(nrow, &rows);
10207: nrow = 0;
10208: for (r = 0; r < rEnd-rStart; ++r) if (a[r] == 0.0) rows[nrow++] = r+rStart;
10209: VecRestoreArrayRead(diag, &a);
10210: VecDestroy(&diag);
10211: ISCreateGeneral(PetscObjectComm((PetscObject) mat), nrow, rows, PETSC_OWN_POINTER, is);
10212: } else {
10213: (*mat->ops->findzerodiagonals)(mat, is);
10214: }
10215: return(0);
10216: }
10218: /*@
10219: MatFindOffBlockDiagonalEntries - Finds all the rows of a matrix that have entries outside of the main diagonal block (defined by the matrix block size)
10221: Collective on Mat
10223: Input Parameter:
10224: . mat - the matrix
10226: Output Parameter:
10227: . is - contains the list of rows with off block diagonal entries
10229: Level: developer
10231: Concepts: matrix-vector product
10233: .seealso: MatMultTranspose(), MatMultAdd(), MatMultTransposeAdd()
10234: @*/
10235: PetscErrorCode MatFindOffBlockDiagonalEntries(Mat mat,IS *is)
10236: {
10242: if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
10243: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
10245: if (!mat->ops->findoffblockdiagonalentries) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"This matrix type does not have a find off block diagonal entries defined");
10246: (*mat->ops->findoffblockdiagonalentries)(mat,is);
10247: return(0);
10248: }
10250: /*@C
10251: MatInvertBlockDiagonal - Inverts the block diagonal entries.
10253: Collective on Mat
10255: Input Parameters:
10256: . mat - the matrix
10258: Output Parameters:
10259: . values - the block inverses in column major order (FORTRAN-like)
10261: Note:
10262: This routine is not available from Fortran.
10264: Level: advanced
10266: .seealso: MatInvertBockDiagonalMat
10267: @*/
10268: PetscErrorCode MatInvertBlockDiagonal(Mat mat,const PetscScalar **values)
10269: {
10274: if (!mat->assembled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
10275: if (mat->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
10276: if (!mat->ops->invertblockdiagonal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not supported");
10277: (*mat->ops->invertblockdiagonal)(mat,values);
10278: return(0);
10279: }
10281: /*@
10282: MatInvertBlockDiagonalMat - set matrix C to be the inverted block diagonal of matrix A
10284: Collective on Mat
10286: Input Parameters:
10287: . A - the matrix
10289: Output Parameters:
10290: . C - matrix with inverted block diagonal of A. This matrix should be created and may have its type set.
10292: Level: advanced
10294: .seealso: MatInvertBockDiagonal()
10295: @*/
10296: PetscErrorCode MatInvertBlockDiagonalMat(Mat A,Mat C)
10297: {
10298: PetscErrorCode ierr;
10299: const PetscScalar *vals;
10300: PetscInt *dnnz;
10301: PetscInt M,N,m,n,rstart,rend,bs,i,j;
10304: MatInvertBlockDiagonal(A,&vals);
10305: MatGetBlockSize(A,&bs);
10306: MatGetSize(A,&M,&N);
10307: MatGetLocalSize(A,&m,&n);
10308: MatSetSizes(C,m,n,M,N);
10309: MatSetBlockSize(C,bs);
10310: PetscMalloc1(m/bs,&dnnz);
10311: for(j = 0; j < m/bs; j++) {
10312: dnnz[j] = 1;
10313: }
10314: MatXAIJSetPreallocation(C,bs,dnnz,NULL,NULL,NULL);
10315: PetscFree(dnnz);
10316: MatGetOwnershipRange(C,&rstart,&rend);
10317: MatSetOption(C,MAT_ROW_ORIENTED,PETSC_FALSE);
10318: for (i = rstart/bs; i < rend/bs; i++) {
10319: MatSetValuesBlocked(C,1,&i,1,&i,&vals[(i-rstart/bs)*bs*bs],INSERT_VALUES);
10320: }
10321: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
10322: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
10323: MatSetOption(C,MAT_ROW_ORIENTED,PETSC_TRUE);
10324: return(0);
10325: }
10327: /*@C
10328: MatTransposeColoringDestroy - Destroys a coloring context for matrix product C=A*B^T that was created
10329: via MatTransposeColoringCreate().
10331: Collective on MatTransposeColoring
10333: Input Parameter:
10334: . c - coloring context
10336: Level: intermediate
10338: .seealso: MatTransposeColoringCreate()
10339: @*/
10340: PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring *c)
10341: {
10342: PetscErrorCode ierr;
10343: MatTransposeColoring matcolor=*c;
10346: if (!matcolor) return(0);
10347: if (--((PetscObject)matcolor)->refct > 0) {matcolor = 0; return(0);}
10349: PetscFree3(matcolor->ncolumns,matcolor->nrows,matcolor->colorforrow);
10350: PetscFree(matcolor->rows);
10351: PetscFree(matcolor->den2sp);
10352: PetscFree(matcolor->colorforcol);
10353: PetscFree(matcolor->columns);
10354: if (matcolor->brows>0) {
10355: PetscFree(matcolor->lstart);
10356: }
10357: PetscHeaderDestroy(c);
10358: return(0);
10359: }
10361: /*@C
10362: MatTransColoringApplySpToDen - Given a symbolic matrix product C=A*B^T for which
10363: a MatTransposeColoring context has been created, computes a dense B^T by Apply
10364: MatTransposeColoring to sparse B.
10366: Collective on MatTransposeColoring
10368: Input Parameters:
10369: + B - sparse matrix B
10370: . Btdense - symbolic dense matrix B^T
10371: - coloring - coloring context created with MatTransposeColoringCreate()
10373: Output Parameter:
10374: . Btdense - dense matrix B^T
10376: Level: advanced
10378: Notes: These are used internally for some implementations of MatRARt()
10380: .seealso: MatTransposeColoringCreate(), MatTransposeColoringDestroy(), MatTransColoringApplyDenToSp()
10382: .keywords: coloring
10383: @*/
10384: PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring coloring,Mat B,Mat Btdense)
10385: {
10393: if (!B->ops->transcoloringapplysptoden) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not supported for this matrix type %s",((PetscObject)B)->type_name);
10394: (B->ops->transcoloringapplysptoden)(coloring,B,Btdense);
10395: return(0);
10396: }
10398: /*@C
10399: MatTransColoringApplyDenToSp - Given a symbolic matrix product Csp=A*B^T for which
10400: a MatTransposeColoring context has been created and a dense matrix Cden=A*Btdense
10401: in which Btdens is obtained from MatTransColoringApplySpToDen(), recover sparse matrix
10402: Csp from Cden.
10404: Collective on MatTransposeColoring
10406: Input Parameters:
10407: + coloring - coloring context created with MatTransposeColoringCreate()
10408: - Cden - matrix product of a sparse matrix and a dense matrix Btdense
10410: Output Parameter:
10411: . Csp - sparse matrix
10413: Level: advanced
10415: Notes: These are used internally for some implementations of MatRARt()
10417: .seealso: MatTransposeColoringCreate(), MatTransposeColoringDestroy(), MatTransColoringApplySpToDen()
10419: .keywords: coloring
10420: @*/
10421: PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring matcoloring,Mat Cden,Mat Csp)
10422: {
10430: if (!Csp->ops->transcoloringapplydentosp) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not supported for this matrix type %s",((PetscObject)Csp)->type_name);
10431: (Csp->ops->transcoloringapplydentosp)(matcoloring,Cden,Csp);
10432: return(0);
10433: }
10435: /*@C
10436: MatTransposeColoringCreate - Creates a matrix coloring context for matrix product C=A*B^T.
10438: Collective on Mat
10440: Input Parameters:
10441: + mat - the matrix product C
10442: - iscoloring - the coloring of the matrix; usually obtained with MatColoringCreate() or DMCreateColoring()
10444: Output Parameter:
10445: . color - the new coloring context
10447: Level: intermediate
10449: .seealso: MatTransposeColoringDestroy(), MatTransColoringApplySpToDen(),
10450: MatTransColoringApplyDenToSp()
10451: @*/
10452: PetscErrorCode MatTransposeColoringCreate(Mat mat,ISColoring iscoloring,MatTransposeColoring *color)
10453: {
10454: MatTransposeColoring c;
10455: MPI_Comm comm;
10456: PetscErrorCode ierr;
10459: PetscLogEventBegin(MAT_TransposeColoringCreate,mat,0,0,0);
10460: PetscObjectGetComm((PetscObject)mat,&comm);
10461: PetscHeaderCreate(c,MAT_TRANSPOSECOLORING_CLASSID,"MatTransposeColoring","Matrix product C=A*B^T via coloring","Mat",comm,MatTransposeColoringDestroy,NULL);
10463: c->ctype = iscoloring->ctype;
10464: if (mat->ops->transposecoloringcreate) {
10465: (*mat->ops->transposecoloringcreate)(mat,iscoloring,c);
10466: } else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Code not yet written for this matrix type");
10468: *color = c;
10469: PetscLogEventEnd(MAT_TransposeColoringCreate,mat,0,0,0);
10470: return(0);
10471: }
10473: /*@
10474: MatGetNonzeroState - Returns a 64 bit integer representing the current state of nonzeros in the matrix. If the
10475: matrix has had no new nonzero locations added to the matrix since the previous call then the value will be the
10476: same, otherwise it will be larger
10478: Not Collective
10480: Input Parameter:
10481: . A - the matrix
10483: Output Parameter:
10484: . state - the current state
10486: Notes: You can only compare states from two different calls to the SAME matrix, you cannot compare calls between
10487: different matrices
10489: Level: intermediate
10491: @*/
10492: PetscErrorCode MatGetNonzeroState(Mat mat,PetscObjectState *state)
10493: {
10496: *state = mat->nonzerostate;
10497: return(0);
10498: }
10500: /*@
10501: MatCreateMPIMatConcatenateSeqMat - Creates a single large PETSc matrix by concatenating sequential
10502: matrices from each processor
10504: Collective on MPI_Comm
10506: Input Parameters:
10507: + comm - the communicators the parallel matrix will live on
10508: . seqmat - the input sequential matrices
10509: . n - number of local columns (or PETSC_DECIDE)
10510: - reuse - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
10512: Output Parameter:
10513: . mpimat - the parallel matrix generated
10515: Level: advanced
10517: Notes: The number of columns of the matrix in EACH processor MUST be the same.
10519: @*/
10520: PetscErrorCode MatCreateMPIMatConcatenateSeqMat(MPI_Comm comm,Mat seqmat,PetscInt n,MatReuse reuse,Mat *mpimat)
10521: {
10525: if (!seqmat->ops->creatempimatconcatenateseqmat) SETERRQ1(PetscObjectComm((PetscObject)seqmat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)seqmat)->type_name);
10526: if (reuse == MAT_REUSE_MATRIX && seqmat == *mpimat) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"MAT_REUSE_MATRIX means reuse the matrix passed in as the final argument, not the original matrix");
10528: PetscLogEventBegin(MAT_Merge,seqmat,0,0,0);
10529: (*seqmat->ops->creatempimatconcatenateseqmat)(comm,seqmat,n,reuse,mpimat);
10530: PetscLogEventEnd(MAT_Merge,seqmat,0,0,0);
10531: return(0);
10532: }
10534: /*@
10535: MatSubdomainsCreateCoalesce - Creates index subdomains by coalescing adjacent
10536: ranks' ownership ranges.
10538: Collective on A
10540: Input Parameters:
10541: + A - the matrix to create subdomains from
10542: - N - requested number of subdomains
10545: Output Parameters:
10546: + n - number of subdomains resulting on this rank
10547: - iss - IS list with indices of subdomains on this rank
10549: Level: advanced
10551: Notes: number of subdomains must be smaller than the communicator size
10552: @*/
10553: PetscErrorCode MatSubdomainsCreateCoalesce(Mat A,PetscInt N,PetscInt *n,IS *iss[])
10554: {
10555: MPI_Comm comm,subcomm;
10556: PetscMPIInt size,rank,color;
10557: PetscInt rstart,rend,k;
10558: PetscErrorCode ierr;
10561: PetscObjectGetComm((PetscObject)A,&comm);
10562: MPI_Comm_size(comm,&size);
10563: MPI_Comm_rank(comm,&rank);
10564: if (N < 1 || N >= (PetscInt)size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of subdomains must be > 0 and < %D, got N = %D",size,N);
10565: *n = 1;
10566: k = ((PetscInt)size)/N + ((PetscInt)size%N>0); /* There are up to k ranks to a color */
10567: color = rank/k;
10568: MPI_Comm_split(comm,color,rank,&subcomm);
10569: PetscMalloc1(1,iss);
10570: MatGetOwnershipRange(A,&rstart,&rend);
10571: ISCreateStride(subcomm,rend-rstart,rstart,1,iss[0]);
10572: MPI_Comm_free(&subcomm);
10573: return(0);
10574: }
10576: /*@
10577: MatGalerkin - Constructs the coarse grid problem via Galerkin projection.
10579: If the interpolation and restriction operators are the same, uses MatPtAP.
10580: If they are not the same, use MatMatMatMult.
10582: Once the coarse grid problem is constructed, correct for interpolation operators
10583: that are not of full rank, which can legitimately happen in the case of non-nested
10584: geometric multigrid.
10586: Input Parameters:
10587: + restrct - restriction operator
10588: . dA - fine grid matrix
10589: . interpolate - interpolation operator
10590: . reuse - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
10591: - fill - expected fill, use PETSC_DEFAULT if you do not have a good estimate
10593: Output Parameters:
10594: . A - the Galerkin coarse matrix
10596: Options Database Key:
10597: . -pc_mg_galerkin <both,pmat,mat,none>
10599: Level: developer
10601: .keywords: MG, multigrid, Galerkin
10603: .seealso: MatPtAP(), MatMatMatMult()
10604: @*/
10605: PetscErrorCode MatGalerkin(Mat restrct, Mat dA, Mat interpolate, MatReuse reuse, PetscReal fill, Mat *A)
10606: {
10608: IS zerorows;
10609: Vec diag;
10612: if (reuse == MAT_INPLACE_MATRIX) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Inplace product not supported");
10613: /* Construct the coarse grid matrix */
10614: if (interpolate == restrct) {
10615: MatPtAP(dA,interpolate,reuse,fill,A);
10616: } else {
10617: MatMatMatMult(restrct,dA,interpolate,reuse,fill,A);
10618: }
10620: /* If the interpolation matrix is not of full rank, A will have zero rows.
10621: This can legitimately happen in the case of non-nested geometric multigrid.
10622: In that event, we set the rows of the matrix to the rows of the identity,
10623: ignoring the equations (as the RHS will also be zero). */
10625: MatFindZeroRows(*A, &zerorows);
10627: if (zerorows != NULL) { /* if there are any zero rows */
10628: MatCreateVecs(*A, &diag, NULL);
10629: MatGetDiagonal(*A, diag);
10630: VecISSet(diag, zerorows, 1.0);
10631: MatDiagonalSet(*A, diag, INSERT_VALUES);
10632: VecDestroy(&diag);
10633: ISDestroy(&zerorows);
10634: }
10635: return(0);
10636: }
10638: /*@C
10639: MatSetOperation - Allows user to set a matrix operation for any matrix type
10641: Logically Collective on Mat
10643: Input Parameters:
10644: + mat - the matrix
10645: . op - the name of the operation
10646: - f - the function that provides the operation
10648: Level: developer
10650: Usage:
10651: $ extern PetscErrorCode usermult(Mat,Vec,Vec);
10652: $ MatCreateXXX(comm,...&A);
10653: $ MatSetOperation(A,MATOP_MULT,(void(*)(void))usermult);
10655: Notes:
10656: See the file include/petscmat.h for a complete list of matrix
10657: operations, which all have the form MATOP_<OPERATION>, where
10658: <OPERATION> is the name (in all capital letters) of the
10659: user interface routine (e.g., MatMult() -> MATOP_MULT).
10661: All user-provided functions (except for MATOP_DESTROY) should have the same calling
10662: sequence as the usual matrix interface routines, since they
10663: are intended to be accessed via the usual matrix interface
10664: routines, e.g.,
10665: $ MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)
10667: In particular each function MUST return an error code of 0 on success and
10668: nonzero on failure.
10670: This routine is distinct from MatShellSetOperation() in that it can be called on any matrix type.
10672: .keywords: matrix, set, operation
10674: .seealso: MatGetOperation(), MatCreateShell(), MatShellSetContext(), MatShellSetOperation()
10675: @*/
10676: PetscErrorCode MatSetOperation(Mat mat,MatOperation op,void (*f)(void))
10677: {
10680: if (op == MATOP_VIEW && !mat->ops->viewnative && f != (void (*)(void))(mat->ops->view)) {
10681: mat->ops->viewnative = mat->ops->view;
10682: }
10683: (((void(**)(void))mat->ops)[op]) = f;
10684: return(0);
10685: }
10687: /*@C
10688: MatGetOperation - Gets a matrix operation for any matrix type.
10690: Not Collective
10692: Input Parameters:
10693: + mat - the matrix
10694: - op - the name of the operation
10696: Output Parameter:
10697: . f - the function that provides the operation
10699: Level: developer
10701: Usage:
10702: $ PetscErrorCode (*usermult)(Mat,Vec,Vec);
10703: $ MatGetOperation(A,MATOP_MULT,(void(**)(void))&usermult);
10705: Notes:
10706: See the file include/petscmat.h for a complete list of matrix
10707: operations, which all have the form MATOP_<OPERATION>, where
10708: <OPERATION> is the name (in all capital letters) of the
10709: user interface routine (e.g., MatMult() -> MATOP_MULT).
10711: This routine is distinct from MatShellGetOperation() in that it can be called on any matrix type.
10713: .keywords: matrix, get, operation
10715: .seealso: MatSetOperation(), MatCreateShell(), MatShellGetContext(), MatShellGetOperation()
10716: @*/
10717: PetscErrorCode MatGetOperation(Mat mat,MatOperation op,void(**f)(void))
10718: {
10721: *f = (((void (**)(void))mat->ops)[op]);
10722: return(0);
10723: }
10725: /*@
10726: MatHasOperation - Determines whether the given matrix supports the particular
10727: operation.
10729: Not Collective
10731: Input Parameters:
10732: + mat - the matrix
10733: - op - the operation, for example, MATOP_GET_DIAGONAL
10735: Output Parameter:
10736: . has - either PETSC_TRUE or PETSC_FALSE
10738: Level: advanced
10740: Notes:
10741: See the file include/petscmat.h for a complete list of matrix
10742: operations, which all have the form MATOP_<OPERATION>, where
10743: <OPERATION> is the name (in all capital letters) of the
10744: user-level routine. E.g., MatNorm() -> MATOP_NORM.
10746: .keywords: matrix, has, operation
10748: .seealso: MatCreateShell()
10749: @*/
10750: PetscErrorCode MatHasOperation(Mat mat,MatOperation op,PetscBool *has)
10751: {
10758: if (mat->ops->hasoperation) {
10759: (*mat->ops->hasoperation)(mat,op,has);
10760: } else {
10761: if (((void**)mat->ops)[op]) *has = PETSC_TRUE;
10762: else {
10763: *has = PETSC_FALSE;
10764: if (op == MATOP_CREATE_SUBMATRIX) {
10765: PetscMPIInt size;
10767: MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);
10768: if (size == 1) {
10769: MatHasOperation(mat,MATOP_CREATE_SUBMATRICES,has);
10770: }
10771: }
10772: }
10773: }
10774: return(0);
10775: }