Actual source code: matrix.c

  1: /*
  2:    This is where the abstract matrix operations are defined
  3: */

  5: #include <petsc/private/matimpl.h>
  6: #include <petsc/private/isimpl.h>
  7: #include <petsc/private/vecimpl.h>

  9: /* Logging support */
 10: PetscClassId MAT_CLASSID;
 11: PetscClassId MAT_COLORING_CLASSID;
 12: PetscClassId MAT_FDCOLORING_CLASSID;
 13: PetscClassId MAT_TRANSPOSECOLORING_CLASSID;

 15: PetscLogEvent MAT_Mult, MAT_Mults, MAT_MultConstrained, MAT_MultAdd, MAT_MultTranspose;
 16: PetscLogEvent MAT_MultTransposeConstrained, MAT_MultTransposeAdd, MAT_Solve, MAT_Solves, MAT_SolveAdd, MAT_SolveTranspose, MAT_MatSolve,MAT_MatTrSolve;
 17: PetscLogEvent MAT_SolveTransposeAdd, MAT_SOR, MAT_ForwardSolve, MAT_BackwardSolve, MAT_LUFactor, MAT_LUFactorSymbolic;
 18: PetscLogEvent MAT_LUFactorNumeric, MAT_CholeskyFactor, MAT_CholeskyFactorSymbolic, MAT_CholeskyFactorNumeric, MAT_ILUFactor;
 19: PetscLogEvent MAT_ILUFactorSymbolic, MAT_ICCFactorSymbolic, MAT_Copy, MAT_Convert, MAT_Scale, MAT_AssemblyBegin;
 20: PetscLogEvent MAT_QRFactorNumeric, MAT_QRFactorSymbolic, MAT_QRFactor;
 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_PartitioningND, 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_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_CUSPARSECopyFromGPU, MAT_CUSPARSEGenerateTranspose, MAT_CUSPARSESolveAnalysis;
 36: PetscLogEvent MAT_PreallCOO, MAT_SetVCOO;
 37: PetscLogEvent MAT_SetValuesBatch;
 38: PetscLogEvent MAT_ViennaCLCopyToGPU;
 39: PetscLogEvent MAT_DenseCopyToGPU, MAT_DenseCopyFromGPU;
 40: PetscLogEvent MAT_Merge,MAT_Residual,MAT_SetRandom;
 41: PetscLogEvent MAT_FactorFactS,MAT_FactorInvS;
 42: PetscLogEvent MATCOLORING_Apply,MATCOLORING_Comm,MATCOLORING_Local,MATCOLORING_ISCreate,MATCOLORING_SetUp,MATCOLORING_Weights;

 44: const char *const MatFactorTypes[] = {"NONE","LU","CHOLESKY","ILU","ICC","ILUDT","QR","MatFactorType","MAT_FACTOR_",NULL};

 46: /*@
 47:    MatSetRandom - Sets all components of a matrix to random numbers. For sparse matrices that have been preallocated but not been assembled it randomly selects appropriate locations,
 48:                   for sparse matrices that already have locations it fills the locations with random numbers

 50:    Logically Collective on Mat

 52:    Input Parameters:
 53: +  x  - the matrix
 54: -  rctx - the random number context, formed by PetscRandomCreate(), or NULL and
 55:           it will create one internally.

 57:    Output Parameter:
 58: .  x  - the matrix

 60:    Example of Usage:
 61: .vb
 62:      PetscRandomCreate(PETSC_COMM_WORLD,&rctx);
 63:      MatSetRandom(x,rctx);
 64:      PetscRandomDestroy(rctx);
 65: .ve

 67:    Level: intermediate


 70: .seealso: MatZeroEntries(), MatSetValues(), PetscRandomCreate(), PetscRandomDestroy()
 71: @*/
 72: PetscErrorCode MatSetRandom(Mat x,PetscRandom rctx)
 73: {
 75:   PetscRandom    randObj = NULL;


 82:   if (!x->ops->setrandom) SETERRQ1(PetscObjectComm((PetscObject)x),PETSC_ERR_SUP,"Mat type %s",((PetscObject)x)->type_name);

 84:   if (!rctx) {
 85:     MPI_Comm comm;
 86:     PetscObjectGetComm((PetscObject)x,&comm);
 87:     PetscRandomCreate(comm,&randObj);
 88:     PetscRandomSetFromOptions(randObj);
 89:     rctx = randObj;
 90:   }

 92:   PetscLogEventBegin(MAT_SetRandom,x,rctx,0,0);
 93:   (*x->ops->setrandom)(x,rctx);
 94:   PetscLogEventEnd(MAT_SetRandom,x,rctx,0,0);

 96:   MatAssemblyBegin(x, MAT_FINAL_ASSEMBLY);
 97:   MatAssemblyEnd(x, MAT_FINAL_ASSEMBLY);
 98:   PetscRandomDestroy(&randObj);
 99:   return(0);
100: }

102: /*@
103:    MatFactorGetErrorZeroPivot - returns the pivot value that was determined to be zero and the row it occurred in

105:    Logically Collective on Mat

107:    Input Parameters:
108: .  mat - the factored matrix

110:    Output Parameter:
111: +  pivot - the pivot value computed
112: -  row - the row that the zero pivot occurred. Note that this row must be interpreted carefully due to row reorderings and which processes
113:          the share the matrix

115:    Level: advanced

117:    Notes:
118:     This routine does not work for factorizations done with external packages.

120:     This routine should only be called if MatGetFactorError() returns a value of MAT_FACTOR_NUMERIC_ZEROPIVOT

122:     This can be called on non-factored matrices that come from, for example, matrices used in SOR.

124: .seealso: MatZeroEntries(), MatFactor(), MatGetFactor(), MatLUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatFactorClearError(), MatFactorGetErrorZeroPivot()
125: @*/
126: PetscErrorCode MatFactorGetErrorZeroPivot(Mat mat,PetscReal *pivot,PetscInt *row)
127: {
130:   *pivot = mat->factorerror_zeropivot_value;
131:   *row   = mat->factorerror_zeropivot_row;
132:   return(0);
133: }

135: /*@
136:    MatFactorGetError - gets the error code from a factorization

138:    Logically Collective on Mat

140:    Input Parameters:
141: .  mat - the factored matrix

143:    Output Parameter:
144: .  err  - the error code

146:    Level: advanced

148:    Notes:
149:     This can be called on non-factored matrices that come from, for example, matrices used in SOR.

151: .seealso: MatZeroEntries(), MatFactor(), MatGetFactor(), MatLUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatFactorClearError(), MatFactorGetErrorZeroPivot()
152: @*/
153: PetscErrorCode MatFactorGetError(Mat mat,MatFactorError *err)
154: {
157:   *err = mat->factorerrortype;
158:   return(0);
159: }

161: /*@
162:    MatFactorClearError - clears the error code in a factorization

164:    Logically Collective on Mat

166:    Input Parameter:
167: .  mat - the factored matrix

169:    Level: developer

171:    Notes:
172:     This can be called on non-factored matrices that come from, for example, matrices used in SOR.

174: .seealso: MatZeroEntries(), MatFactor(), MatGetFactor(), MatLUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatFactorGetError(), MatFactorGetErrorZeroPivot()
175: @*/
176: PetscErrorCode MatFactorClearError(Mat mat)
177: {
180:   mat->factorerrortype             = MAT_FACTOR_NOERROR;
181:   mat->factorerror_zeropivot_value = 0.0;
182:   mat->factorerror_zeropivot_row   = 0;
183:   return(0);
184: }

186: PETSC_INTERN PetscErrorCode MatFindNonzeroRowsOrCols_Basic(Mat mat,PetscBool cols,PetscReal tol,IS *nonzero)
187: {
188:   PetscErrorCode    ierr;
189:   Vec               r,l;
190:   const PetscScalar *al;
191:   PetscInt          i,nz,gnz,N,n;

194:   MatCreateVecs(mat,&r,&l);
195:   if (!cols) { /* nonzero rows */
196:     MatGetSize(mat,&N,NULL);
197:     MatGetLocalSize(mat,&n,NULL);
198:     VecSet(l,0.0);
199:     VecSetRandom(r,NULL);
200:     MatMult(mat,r,l);
201:     VecGetArrayRead(l,&al);
202:   } else { /* nonzero columns */
203:     MatGetSize(mat,NULL,&N);
204:     MatGetLocalSize(mat,NULL,&n);
205:     VecSet(r,0.0);
206:     VecSetRandom(l,NULL);
207:     MatMultTranspose(mat,l,r);
208:     VecGetArrayRead(r,&al);
209:   }
210:   if (tol <= 0.0) { for (i=0,nz=0;i<n;i++) if (al[i] != 0.0) nz++; }
211:   else { for (i=0,nz=0;i<n;i++) if (PetscAbsScalar(al[i]) > tol) nz++; }
212:   MPIU_Allreduce(&nz,&gnz,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)mat));
213:   if (gnz != N) {
214:     PetscInt *nzr;
215:     PetscMalloc1(nz,&nzr);
216:     if (nz) {
217:       if (tol < 0) { for (i=0,nz=0;i<n;i++) if (al[i] != 0.0) nzr[nz++] = i; }
218:       else { for (i=0,nz=0;i<n;i++) if (PetscAbsScalar(al[i]) > tol) nzr[nz++] = i; }
219:     }
220:     ISCreateGeneral(PetscObjectComm((PetscObject)mat),nz,nzr,PETSC_OWN_POINTER,nonzero);
221:   } else *nonzero = NULL;
222:   if (!cols) { /* nonzero rows */
223:     VecRestoreArrayRead(l,&al);
224:   } else {
225:     VecRestoreArrayRead(r,&al);
226:   }
227:   VecDestroy(&l);
228:   VecDestroy(&r);
229:   return(0);
230: }

232: /*@
233:       MatFindNonzeroRows - Locate all rows that are not completely zero in the matrix

235:   Input Parameter:
236: .    A  - the matrix

238:   Output Parameter:
239: .    keptrows - the rows that are not completely zero

241:   Notes:
242:     keptrows is set to NULL if all rows are nonzero.

244:   Level: intermediate

246:  @*/
247: PetscErrorCode MatFindNonzeroRows(Mat mat,IS *keptrows)
248: {

255:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
256:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
257:   if (!mat->ops->findnonzerorows) {
258:     MatFindNonzeroRowsOrCols_Basic(mat,PETSC_FALSE,0.0,keptrows);
259:   } else {
260:     (*mat->ops->findnonzerorows)(mat,keptrows);
261:   }
262:   return(0);
263: }

265: /*@
266:       MatFindZeroRows - Locate all rows that are completely zero in the matrix

268:   Input Parameter:
269: .    A  - the matrix

271:   Output Parameter:
272: .    zerorows - the rows that are completely zero

274:   Notes:
275:     zerorows is set to NULL if no rows are zero.

277:   Level: intermediate

279:  @*/
280: PetscErrorCode MatFindZeroRows(Mat mat,IS *zerorows)
281: {
283:   IS keptrows;
284:   PetscInt m, n;


289:   MatFindNonzeroRows(mat, &keptrows);
290:   /* MatFindNonzeroRows sets keptrows to NULL if there are no zero rows.
291:      In keeping with this convention, we set zerorows to NULL if there are no zero
292:      rows. */
293:   if (keptrows == NULL) {
294:     *zerorows = NULL;
295:   } else {
296:     MatGetOwnershipRange(mat,&m,&n);
297:     ISComplement(keptrows,m,n,zerorows);
298:     ISDestroy(&keptrows);
299:   }
300:   return(0);
301: }

303: /*@
304:    MatGetDiagonalBlock - Returns the part of the matrix associated with the on-process coupling

306:    Not Collective

308:    Input Parameters:
309: .   A - the matrix

311:    Output Parameters:
312: .   a - the diagonal part (which is a SEQUENTIAL matrix)

314:    Notes:
315:     see the manual page for MatCreateAIJ() for more information on the "diagonal part" of the matrix.
316:           Use caution, as the reference count on the returned matrix is not incremented and it is used as
317:           part of the containing MPI Mat's normal operation.

319:    Level: advanced

321: @*/
322: PetscErrorCode MatGetDiagonalBlock(Mat A,Mat *a)
323: {

330:   if (A->factortype) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
331:   if (!A->ops->getdiagonalblock) {
332:     PetscMPIInt size;
333:     MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);
334:     if (size == 1) {
335:       *a = A;
336:       return(0);
337:     } else SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Not coded for matrix type %s",((PetscObject)A)->type_name);
338:   }
339:   (*A->ops->getdiagonalblock)(A,a);
340:   return(0);
341: }

343: /*@
344:    MatGetTrace - Gets the trace of a matrix. The sum of the diagonal entries.

346:    Collective on Mat

348:    Input Parameters:
349: .  mat - the matrix

351:    Output Parameter:
352: .   trace - the sum of the diagonal entries

354:    Level: advanced

356: @*/
357: PetscErrorCode MatGetTrace(Mat mat,PetscScalar *trace)
358: {
360:   Vec            diag;

363:   MatCreateVecs(mat,&diag,NULL);
364:   MatGetDiagonal(mat,diag);
365:   VecSum(diag,trace);
366:   VecDestroy(&diag);
367:   return(0);
368: }

370: /*@
371:    MatRealPart - Zeros out the imaginary part of the matrix

373:    Logically Collective on Mat

375:    Input Parameters:
376: .  mat - the matrix

378:    Level: advanced


381: .seealso: MatImaginaryPart()
382: @*/
383: PetscErrorCode MatRealPart(Mat mat)
384: {

390:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
391:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
392:   if (!mat->ops->realpart) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
393:   MatCheckPreallocated(mat,1);
394:   (*mat->ops->realpart)(mat);
395:   return(0);
396: }

398: /*@C
399:    MatGetGhosts - Get the global index of all ghost nodes defined by the sparse matrix

401:    Collective on Mat

403:    Input Parameter:
404: .  mat - the matrix

406:    Output Parameters:
407: +   nghosts - number of ghosts (note for BAIJ matrices there is one ghost for each block)
408: -   ghosts - the global indices of the ghost points

410:    Notes:
411:     the nghosts and ghosts are suitable to pass into VecCreateGhost()

413:    Level: advanced

415: @*/
416: PetscErrorCode MatGetGhosts(Mat mat,PetscInt *nghosts,const PetscInt *ghosts[])
417: {

423:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
424:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
425:   if (!mat->ops->getghosts) {
426:     if (nghosts) *nghosts = 0;
427:     if (ghosts) *ghosts = NULL;
428:   } else {
429:     (*mat->ops->getghosts)(mat,nghosts,ghosts);
430:   }
431:   return(0);
432: }


435: /*@
436:    MatImaginaryPart - Moves the imaginary part of the matrix to the real part and zeros the imaginary part

438:    Logically Collective on Mat

440:    Input Parameters:
441: .  mat - the matrix

443:    Level: advanced


446: .seealso: MatRealPart()
447: @*/
448: PetscErrorCode MatImaginaryPart(Mat mat)
449: {

455:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
456:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
457:   if (!mat->ops->imaginarypart) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
458:   MatCheckPreallocated(mat,1);
459:   (*mat->ops->imaginarypart)(mat);
460:   return(0);
461: }

463: /*@
464:    MatMissingDiagonal - Determine if sparse matrix is missing a diagonal entry (or block entry for BAIJ matrices)

466:    Not Collective

468:    Input Parameter:
469: .  mat - the matrix

471:    Output Parameters:
472: +  missing - is any diagonal missing
473: -  dd - first diagonal entry that is missing (optional) on this process

475:    Level: advanced


478: .seealso: MatRealPart()
479: @*/
480: PetscErrorCode MatMissingDiagonal(Mat mat,PetscBool *missing,PetscInt *dd)
481: {

488:   if (!mat->assembled) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix %s",((PetscObject)mat)->type_name);
489:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
490:   if (!mat->ops->missingdiagonal) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
491:   (*mat->ops->missingdiagonal)(mat,missing,dd);
492:   return(0);
493: }

495: /*@C
496:    MatGetRow - Gets a row of a matrix.  You MUST call MatRestoreRow()
497:    for each row that you get to ensure that your application does
498:    not bleed memory.

500:    Not Collective

502:    Input Parameters:
503: +  mat - the matrix
504: -  row - the row to get

506:    Output Parameters:
507: +  ncols -  if not NULL, the number of nonzeros in the row
508: .  cols - if not NULL, the column numbers
509: -  vals - if not NULL, the values

511:    Notes:
512:    This routine is provided for people who need to have direct access
513:    to the structure of a matrix.  We hope that we provide enough
514:    high-level matrix routines that few users will need it.

516:    MatGetRow() always returns 0-based column indices, regardless of
517:    whether the internal representation is 0-based (default) or 1-based.

519:    For better efficiency, set cols and/or vals to NULL if you do
520:    not wish to extract these quantities.

522:    The user can only examine the values extracted with MatGetRow();
523:    the values cannot be altered.  To change the matrix entries, one
524:    must use MatSetValues().

526:    You can only have one call to MatGetRow() outstanding for a particular
527:    matrix at a time, per processor. MatGetRow() can only obtain rows
528:    associated with the given processor, it cannot get rows from the
529:    other processors; for that we suggest using MatCreateSubMatrices(), then
530:    MatGetRow() on the submatrix. The row index passed to MatGetRow()
531:    is in the global number of rows.

533:    Fortran Notes:
534:    The calling sequence from Fortran is
535: .vb
536:    MatGetRow(matrix,row,ncols,cols,values,ierr)
537:          Mat     matrix (input)
538:          integer row    (input)
539:          integer ncols  (output)
540:          integer cols(maxcols) (output)
541:          double precision (or double complex) values(maxcols) output
542: .ve
543:    where maxcols >= maximum nonzeros in any row of the matrix.


546:    Caution:
547:    Do not try to change the contents of the output arrays (cols and vals).
548:    In some cases, this may corrupt the matrix.

550:    Level: advanced

552: .seealso: MatRestoreRow(), MatSetValues(), MatGetValues(), MatCreateSubMatrices(), MatGetDiagonal()
553: @*/
554: PetscErrorCode MatGetRow(Mat mat,PetscInt row,PetscInt *ncols,const PetscInt *cols[],const PetscScalar *vals[])
555: {
557:   PetscInt       incols;

562:   if (!mat->assembled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
563:   if (mat->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
564:   if (!mat->ops->getrow) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
565:   MatCheckPreallocated(mat,1);
566:   if (row < mat->rmap->rstart || row >= mat->rmap->rend) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Only for local rows, %D not in [%D,%D)",row,mat->rmap->rstart,mat->rmap->rend);
567:   PetscLogEventBegin(MAT_GetRow,mat,0,0,0);
568:   (*mat->ops->getrow)(mat,row,&incols,(PetscInt**)cols,(PetscScalar**)vals);
569:   if (ncols) *ncols = incols;
570:   PetscLogEventEnd(MAT_GetRow,mat,0,0,0);
571:   return(0);
572: }

574: /*@
575:    MatConjugate - replaces the matrix values with their complex conjugates

577:    Logically Collective on Mat

579:    Input Parameters:
580: .  mat - the matrix

582:    Level: advanced

584: .seealso:  VecConjugate()
585: @*/
586: PetscErrorCode MatConjugate(Mat mat)
587: {
588: #if defined(PETSC_USE_COMPLEX)

593:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
594:   if (!mat->ops->conjugate) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Not provided for matrix type %s, send email to petsc-maint@mcs.anl.gov",((PetscObject)mat)->type_name);
595:   (*mat->ops->conjugate)(mat);
596: #else
598: #endif
599:   return(0);
600: }

602: /*@C
603:    MatRestoreRow - Frees any temporary space allocated by MatGetRow().

605:    Not Collective

607:    Input Parameters:
608: +  mat - the matrix
609: .  row - the row to get
610: .  ncols, cols - the number of nonzeros and their columns
611: -  vals - if nonzero the column values

613:    Notes:
614:    This routine should be called after you have finished examining the entries.

616:    This routine zeros out ncols, cols, and vals. This is to prevent accidental
617:    us of the array after it has been restored. If you pass NULL, it will
618:    not zero the pointers.  Use of cols or vals after MatRestoreRow is invalid.

620:    Fortran Notes:
621:    The calling sequence from Fortran is
622: .vb
623:    MatRestoreRow(matrix,row,ncols,cols,values,ierr)
624:       Mat     matrix (input)
625:       integer row    (input)
626:       integer ncols  (output)
627:       integer cols(maxcols) (output)
628:       double precision (or double complex) values(maxcols) output
629: .ve
630:    Where maxcols >= maximum nonzeros in any row of the matrix.

632:    In Fortran MatRestoreRow() MUST be called after MatGetRow()
633:    before another call to MatGetRow() can be made.

635:    Level: advanced

637: .seealso:  MatGetRow()
638: @*/
639: PetscErrorCode MatRestoreRow(Mat mat,PetscInt row,PetscInt *ncols,const PetscInt *cols[],const PetscScalar *vals[])
640: {

646:   if (!mat->assembled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
647:   if (!mat->ops->restorerow) return(0);
648:   (*mat->ops->restorerow)(mat,row,ncols,(PetscInt **)cols,(PetscScalar **)vals);
649:   if (ncols) *ncols = 0;
650:   if (cols)  *cols = NULL;
651:   if (vals)  *vals = NULL;
652:   return(0);
653: }

655: /*@
656:    MatGetRowUpperTriangular - Sets a flag to enable calls to MatGetRow() for matrix in MATSBAIJ format.
657:    You should call MatRestoreRowUpperTriangular() after calling MatGetRow/MatRestoreRow() to disable the flag.

659:    Not Collective

661:    Input Parameters:
662: .  mat - the matrix

664:    Notes:
665:    The flag is to ensure that users are aware of MatGetRow() only provides the upper triangular part of the row for the matrices in MATSBAIJ format.

667:    Level: advanced

669: .seealso: MatRestoreRowUpperTriangular()
670: @*/
671: PetscErrorCode MatGetRowUpperTriangular(Mat mat)
672: {

678:   if (!mat->assembled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
679:   if (mat->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
680:   MatCheckPreallocated(mat,1);
681:   if (!mat->ops->getrowuppertriangular) return(0);
682:   (*mat->ops->getrowuppertriangular)(mat);
683:   return(0);
684: }

686: /*@
687:    MatRestoreRowUpperTriangular - Disable calls to MatGetRow() for matrix in MATSBAIJ format.

689:    Not Collective

691:    Input Parameters:
692: .  mat - the matrix

694:    Notes:
695:    This routine should be called after you have finished MatGetRow/MatRestoreRow().


698:    Level: advanced

700: .seealso:  MatGetRowUpperTriangular()
701: @*/
702: PetscErrorCode MatRestoreRowUpperTriangular(Mat mat)
703: {

709:   if (!mat->assembled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
710:   if (mat->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
711:   MatCheckPreallocated(mat,1);
712:   if (!mat->ops->restorerowuppertriangular) return(0);
713:   (*mat->ops->restorerowuppertriangular)(mat);
714:   return(0);
715: }

717: /*@C
718:    MatSetOptionsPrefix - Sets the prefix used for searching for all
719:    Mat options in the database.

721:    Logically Collective on Mat

723:    Input Parameter:
724: +  A - the Mat context
725: -  prefix - the prefix to prepend to all option names

727:    Notes:
728:    A hyphen (-) must NOT be given at the beginning of the prefix name.
729:    The first character of all runtime options is AUTOMATICALLY the hyphen.

731:    Level: advanced

733: .seealso: MatSetFromOptions()
734: @*/
735: PetscErrorCode MatSetOptionsPrefix(Mat A,const char prefix[])
736: {

741:   PetscObjectSetOptionsPrefix((PetscObject)A,prefix);
742:   return(0);
743: }

745: /*@C
746:    MatAppendOptionsPrefix - Appends to the prefix used for searching for all
747:    Mat options in the database.

749:    Logically Collective on Mat

751:    Input Parameters:
752: +  A - the Mat context
753: -  prefix - the prefix to prepend to all option names

755:    Notes:
756:    A hyphen (-) must NOT be given at the beginning of the prefix name.
757:    The first character of all runtime options is AUTOMATICALLY the hyphen.

759:    Level: advanced

761: .seealso: MatGetOptionsPrefix()
762: @*/
763: PetscErrorCode MatAppendOptionsPrefix(Mat A,const char prefix[])
764: {

769:   PetscObjectAppendOptionsPrefix((PetscObject)A,prefix);
770:   return(0);
771: }

773: /*@C
774:    MatGetOptionsPrefix - Gets the prefix used for searching for all
775:    Mat options in the database.

777:    Not Collective

779:    Input Parameter:
780: .  A - the Mat context

782:    Output Parameter:
783: .  prefix - pointer to the prefix string used

785:    Notes:
786:     On the fortran side, the user should pass in a string 'prefix' of
787:    sufficient length to hold the prefix.

789:    Level: advanced

791: .seealso: MatAppendOptionsPrefix()
792: @*/
793: PetscErrorCode MatGetOptionsPrefix(Mat A,const char *prefix[])
794: {

799:   PetscObjectGetOptionsPrefix((PetscObject)A,prefix);
800:   return(0);
801: }

803: /*@
804:    MatResetPreallocation - Reset mat to use the original nonzero pattern provided by users.

806:    Collective on Mat

808:    Input Parameters:
809: .  A - the Mat context

811:    Notes:
812:    The allocated memory will be shrunk after calling MatAssembly with MAT_FINAL_ASSEMBLY. Users can reset the preallocation to access the original memory.
813:    Currently support MPIAIJ and SEQAIJ.

815:    Level: beginner

817: .seealso: MatSeqAIJSetPreallocation(), MatMPIAIJSetPreallocation(), MatXAIJSetPreallocation()
818: @*/
819: PetscErrorCode MatResetPreallocation(Mat A)
820: {

826:   PetscUseMethod(A,"MatResetPreallocation_C",(Mat),(A));
827:   return(0);
828: }


831: /*@
832:    MatSetUp - Sets up the internal matrix data structures for later use.

834:    Collective on Mat

836:    Input Parameters:
837: .  A - the Mat context

839:    Notes:
840:    If the user has not set preallocation for this matrix then a default preallocation that is likely to be inefficient is used.

842:    If a suitable preallocation routine is used, this function does not need to be called.

844:    See the Performance chapter of the PETSc users manual for how to preallocate matrices

846:    Level: beginner

848: .seealso: MatCreate(), MatDestroy()
849: @*/
850: PetscErrorCode MatSetUp(Mat A)
851: {
852:   PetscMPIInt    size;

857:   if (!((PetscObject)A)->type_name) {
858:     MPI_Comm_size(PetscObjectComm((PetscObject)A), &size);
859:     if (size == 1) {
860:       MatSetType(A, MATSEQAIJ);
861:     } else {
862:       MatSetType(A, MATMPIAIJ);
863:     }
864:   }
865:   if (!A->preallocated && A->ops->setup) {
866:     PetscInfo(A,"Warning not preallocating matrix storage\n");
867:     (*A->ops->setup)(A);
868:   }
869:   PetscLayoutSetUp(A->rmap);
870:   PetscLayoutSetUp(A->cmap);
871:   A->preallocated = PETSC_TRUE;
872:   return(0);
873: }

875: #if defined(PETSC_HAVE_SAWS)
876: #include <petscviewersaws.h>
877: #endif

879: /*@C
880:    MatViewFromOptions - View from Options

882:    Collective on Mat

884:    Input Parameters:
885: +  A - the Mat context
886: .  obj - Optional object
887: -  name - command line option

889:    Level: intermediate
890: .seealso:  Mat, MatView, PetscObjectViewFromOptions(), MatCreate()
891: @*/
892: PetscErrorCode  MatViewFromOptions(Mat A,PetscObject obj,const char name[])
893: {

898:   PetscObjectViewFromOptions((PetscObject)A,obj,name);
899:   return(0);
900: }

902: /*@C
903:    MatView - Visualizes a matrix object.

905:    Collective on Mat

907:    Input Parameters:
908: +  mat - the matrix
909: -  viewer - visualization context

911:   Notes:
912:   The available visualization contexts include
913: +    PETSC_VIEWER_STDOUT_SELF - for sequential matrices
914: .    PETSC_VIEWER_STDOUT_WORLD - for parallel matrices created on PETSC_COMM_WORLD
915: .    PETSC_VIEWER_STDOUT_(comm) - for matrices created on MPI communicator comm
916: -     PETSC_VIEWER_DRAW_WORLD - graphical display of nonzero structure

918:    The user can open alternative visualization contexts with
919: +    PetscViewerASCIIOpen() - Outputs matrix to a specified file
920: .    PetscViewerBinaryOpen() - Outputs matrix in binary to a
921:          specified file; corresponding input uses MatLoad()
922: .    PetscViewerDrawOpen() - Outputs nonzero matrix structure to
923:          an X window display
924: -    PetscViewerSocketOpen() - Outputs matrix to Socket viewer.
925:          Currently only the sequential dense and AIJ
926:          matrix types support the Socket viewer.

928:    The user can call PetscViewerPushFormat() to specify the output
929:    format of ASCII printed objects (when using PETSC_VIEWER_STDOUT_SELF,
930:    PETSC_VIEWER_STDOUT_WORLD and PetscViewerASCIIOpen).  Available formats include
931: +    PETSC_VIEWER_DEFAULT - default, prints matrix contents
932: .    PETSC_VIEWER_ASCII_MATLAB - prints matrix contents in Matlab format
933: .    PETSC_VIEWER_ASCII_DENSE - prints entire matrix including zeros
934: .    PETSC_VIEWER_ASCII_COMMON - prints matrix contents, using a sparse
935:          format common among all matrix types
936: .    PETSC_VIEWER_ASCII_IMPL - prints matrix contents, using an implementation-specific
937:          format (which is in many cases the same as the default)
938: .    PETSC_VIEWER_ASCII_INFO - prints basic information about the matrix
939:          size and structure (not the matrix entries)
940: -    PETSC_VIEWER_ASCII_INFO_DETAIL - prints more detailed information about
941:          the matrix structure

943:    Options Database Keys:
944: +  -mat_view ::ascii_info - Prints info on matrix at conclusion of MatAssemblyEnd()
945: .  -mat_view ::ascii_info_detail - Prints more detailed info
946: .  -mat_view - Prints matrix in ASCII format
947: .  -mat_view ::ascii_matlab - Prints matrix in Matlab format
948: .  -mat_view draw - PetscDraws nonzero structure of matrix, using MatView() and PetscDrawOpenX().
949: .  -display <name> - Sets display name (default is host)
950: .  -draw_pause <sec> - Sets number of seconds to pause after display
951: .  -mat_view socket - Sends matrix to socket, can be accessed from Matlab (see Users-Manual: ch_matlab for details)
952: .  -viewer_socket_machine <machine> -
953: .  -viewer_socket_port <port> -
954: .  -mat_view binary - save matrix to file in binary format
955: -  -viewer_binary_filename <name> -
956:    Level: beginner

958:    Notes:
959:     The ASCII viewers are only recommended for small matrices on at most a moderate number of processes,
960:     the program will seemingly hang and take hours for larger matrices, for larger matrices one should use the binary format.

962:     In the debugger you can do "call MatView(mat,0)" to display the matrix. (The same holds for any PETSc object viewer).

964:     See the manual page for MatLoad() for the exact format of the binary file when the binary
965:       viewer is used.

967:       See share/petsc/matlab/PetscBinaryRead.m for a Matlab code that can read in the binary file when the binary
968:       viewer is used and lib/petsc/bin/PetscBinaryIO.py for loading them into Python.

970:       One can use '-mat_view draw -draw_pause -1' to pause the graphical display of matrix nonzero structure,
971:       and then use the following mouse functions.
972: + left mouse: zoom in
973: . middle mouse: zoom out
974: - right mouse: continue with the simulation

976: .seealso: PetscViewerPushFormat(), PetscViewerASCIIOpen(), PetscViewerDrawOpen(),
977:           PetscViewerSocketOpen(), PetscViewerBinaryOpen(), MatLoad()
978: @*/
979: PetscErrorCode MatView(Mat mat,PetscViewer viewer)
980: {
981:   PetscErrorCode    ierr;
982:   PetscInt          rows,cols,rbs,cbs;
983:   PetscBool         isascii,isstring,issaws;
984:   PetscViewerFormat format;
985:   PetscMPIInt       size;

990:   if (!viewer) {PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)mat),&viewer);}
993:   MatCheckPreallocated(mat,1);

995:   PetscViewerGetFormat(viewer,&format);
996:   MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);
997:   if (size == 1 && format == PETSC_VIEWER_LOAD_BALANCE) return(0);

999:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
1000:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
1001:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSAWS,&issaws);
1002:   if ((!isascii || (format != PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL)) && mat->factortype) {
1003:     SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"No viewers for factored matrix except ASCII info or info_detail");
1004:   }

1006:   PetscLogEventBegin(MAT_View,mat,viewer,0,0);
1007:   if (isascii) {
1008:     if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ORDER,"Must call MatAssemblyBegin/End() before viewing matrix");
1009:     PetscObjectPrintClassNamePrefixType((PetscObject)mat,viewer);
1010:     if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
1011:       MatNullSpace nullsp,transnullsp;

1013:       PetscViewerASCIIPushTab(viewer);
1014:       MatGetSize(mat,&rows,&cols);
1015:       MatGetBlockSizes(mat,&rbs,&cbs);
1016:       if (rbs != 1 || cbs != 1) {
1017:         if (rbs != cbs) {PetscViewerASCIIPrintf(viewer,"rows=%D, cols=%D, rbs=%D, cbs=%D\n",rows,cols,rbs,cbs);}
1018:         else            {PetscViewerASCIIPrintf(viewer,"rows=%D, cols=%D, bs=%D\n",rows,cols,rbs);}
1019:       } else {
1020:         PetscViewerASCIIPrintf(viewer,"rows=%D, cols=%D\n",rows,cols);
1021:       }
1022:       if (mat->factortype) {
1023:         MatSolverType solver;
1024:         MatFactorGetSolverType(mat,&solver);
1025:         PetscViewerASCIIPrintf(viewer,"package used to perform factorization: %s\n",solver);
1026:       }
1027:       if (mat->ops->getinfo) {
1028:         MatInfo info;
1029:         MatGetInfo(mat,MAT_GLOBAL_SUM,&info);
1030:         PetscViewerASCIIPrintf(viewer,"total: nonzeros=%.f, allocated nonzeros=%.f\n",info.nz_used,info.nz_allocated);
1031:         if (!mat->factortype) {
1032:           PetscViewerASCIIPrintf(viewer,"total number of mallocs used during MatSetValues calls=%D\n",(PetscInt)info.mallocs);
1033:         }
1034:       }
1035:       MatGetNullSpace(mat,&nullsp);
1036:       MatGetTransposeNullSpace(mat,&transnullsp);
1037:       if (nullsp) {PetscViewerASCIIPrintf(viewer,"  has attached null space\n");}
1038:       if (transnullsp && transnullsp != nullsp) {PetscViewerASCIIPrintf(viewer,"  has attached transposed null space\n");}
1039:       MatGetNearNullSpace(mat,&nullsp);
1040:       if (nullsp) {PetscViewerASCIIPrintf(viewer,"  has attached near null space\n");}
1041:       PetscViewerASCIIPushTab(viewer);
1042:       MatProductView(mat,viewer);
1043:       PetscViewerASCIIPopTab(viewer);
1044:     }
1045:   } else if (issaws) {
1046: #if defined(PETSC_HAVE_SAWS)
1047:     PetscMPIInt rank;

1049:     PetscObjectName((PetscObject)mat);
1050:     MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
1051:     if (!((PetscObject)mat)->amsmem && !rank) {
1052:       PetscObjectViewSAWs((PetscObject)mat,viewer);
1053:     }
1054: #endif
1055:   } else if (isstring) {
1056:     const char *type;
1057:     MatGetType(mat,&type);
1058:     PetscViewerStringSPrintf(viewer," MatType: %-7.7s",type);
1059:     if (mat->ops->view) {(*mat->ops->view)(mat,viewer);}
1060:   }
1061:   if ((format == PETSC_VIEWER_NATIVE || format == PETSC_VIEWER_LOAD_BALANCE) && mat->ops->viewnative) {
1062:     PetscViewerASCIIPushTab(viewer);
1063:     (*mat->ops->viewnative)(mat,viewer);
1064:     PetscViewerASCIIPopTab(viewer);
1065:   } else if (mat->ops->view) {
1066:     PetscViewerASCIIPushTab(viewer);
1067:     (*mat->ops->view)(mat,viewer);
1068:     PetscViewerASCIIPopTab(viewer);
1069:   }
1070:   if (isascii) {
1071:     PetscViewerGetFormat(viewer,&format);
1072:     if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
1073:       PetscViewerASCIIPopTab(viewer);
1074:     }
1075:   }
1076:   PetscLogEventEnd(MAT_View,mat,viewer,0,0);
1077:   return(0);
1078: }

1080: #if defined(PETSC_USE_DEBUG)
1081: #include <../src/sys/totalview/tv_data_display.h>
1082: PETSC_UNUSED static int TV_display_type(const struct _p_Mat *mat)
1083: {
1084:   TV_add_row("Local rows", "int", &mat->rmap->n);
1085:   TV_add_row("Local columns", "int", &mat->cmap->n);
1086:   TV_add_row("Global rows", "int", &mat->rmap->N);
1087:   TV_add_row("Global columns", "int", &mat->cmap->N);
1088:   TV_add_row("Typename", TV_ascii_string_type, ((PetscObject)mat)->type_name);
1089:   return TV_format_OK;
1090: }
1091: #endif

1093: /*@C
1094:    MatLoad - Loads a matrix that has been stored in binary/HDF5 format
1095:    with MatView().  The matrix format is determined from the options database.
1096:    Generates a parallel MPI matrix if the communicator has more than one
1097:    processor.  The default matrix type is AIJ.

1099:    Collective on PetscViewer

1101:    Input Parameters:
1102: +  mat - the newly loaded matrix, this needs to have been created with MatCreate()
1103:             or some related function before a call to MatLoad()
1104: -  viewer - binary/HDF5 file viewer

1106:    Options Database Keys:
1107:    Used with block matrix formats (MATSEQBAIJ,  ...) to specify
1108:    block size
1109: .    -matload_block_size <bs>

1111:    Level: beginner

1113:    Notes:
1114:    If the Mat type has not yet been given then MATAIJ is used, call MatSetFromOptions() on the
1115:    Mat before calling this routine if you wish to set it from the options database.

1117:    MatLoad() automatically loads into the options database any options
1118:    given in the file filename.info where filename is the name of the file
1119:    that was passed to the PetscViewerBinaryOpen(). The options in the info
1120:    file will be ignored if you use the -viewer_binary_skip_info option.

1122:    If the type or size of mat is not set before a call to MatLoad, PETSc
1123:    sets the default matrix type AIJ and sets the local and global sizes.
1124:    If type and/or size is already set, then the same are used.

1126:    In parallel, each processor can load a subset of rows (or the
1127:    entire matrix).  This routine is especially useful when a large
1128:    matrix is stored on disk and only part of it is desired on each
1129:    processor.  For example, a parallel solver may access only some of
1130:    the rows from each processor.  The algorithm used here reads
1131:    relatively small blocks of data rather than reading the entire
1132:    matrix and then subsetting it.

1134:    Viewer's PetscViewerType must be either PETSCVIEWERBINARY or PETSCVIEWERHDF5.
1135:    Such viewer can be created using PetscViewerBinaryOpen()/PetscViewerHDF5Open(),
1136:    or the sequence like
1137: $    PetscViewer v;
1138: $    PetscViewerCreate(PETSC_COMM_WORLD,&v);
1139: $    PetscViewerSetType(v,PETSCVIEWERBINARY);
1140: $    PetscViewerSetFromOptions(v);
1141: $    PetscViewerFileSetMode(v,FILE_MODE_READ);
1142: $    PetscViewerFileSetName(v,"datafile");
1143:    The optional PetscViewerSetFromOptions() call allows to override PetscViewerSetType() using option
1144: $ -viewer_type {binary,hdf5}

1146:    See the example src/ksp/ksp/tutorials/ex27.c with the first approach,
1147:    and src/mat/tutorials/ex10.c with the second approach.

1149:    Notes about the PETSc binary format:
1150:    In case of PETSCVIEWERBINARY, a native PETSc binary format is used. Each of the blocks
1151:    is read onto rank 0 and then shipped to its destination rank, one after another.
1152:    Multiple objects, both matrices and vectors, can be stored within the same file.
1153:    Their PetscObject name is ignored; they are loaded in the order of their storage.

1155:    Most users should not need to know the details of the binary storage
1156:    format, since MatLoad() and MatView() completely hide these details.
1157:    But for anyone who's interested, the standard binary matrix storage
1158:    format is

1160: $    PetscInt    MAT_FILE_CLASSID
1161: $    PetscInt    number of rows
1162: $    PetscInt    number of columns
1163: $    PetscInt    total number of nonzeros
1164: $    PetscInt    *number nonzeros in each row
1165: $    PetscInt    *column indices of all nonzeros (starting index is zero)
1166: $    PetscScalar *values of all nonzeros

1168:    PETSc automatically does the byte swapping for
1169: machines that store the bytes reversed, e.g.  DEC alpha, freebsd,
1170: linux, Windows and the paragon; thus if you write your own binary
1171: read/write routines you have to swap the bytes; see PetscBinaryRead()
1172: and PetscBinaryWrite() to see how this may be done.

1174:    Notes about the HDF5 (MATLAB MAT-File Version 7.3) format:
1175:    In case of PETSCVIEWERHDF5, a parallel HDF5 reader is used.
1176:    Each processor's chunk is loaded independently by its owning rank.
1177:    Multiple objects, both matrices and vectors, can be stored within the same file.
1178:    They are looked up by their PetscObject name.

1180:    As the MATLAB MAT-File Version 7.3 format is also a HDF5 flavor, we decided to use
1181:    by default the same structure and naming of the AIJ arrays and column count
1182:    within the HDF5 file. This means that a MAT file saved with -v7.3 flag, e.g.
1183: $    save example.mat A b -v7.3
1184:    can be directly read by this routine (see Reference 1 for details).
1185:    Note that depending on your MATLAB version, this format might be a default,
1186:    otherwise you can set it as default in Preferences.

1188:    Unless -nocompression flag is used to save the file in MATLAB,
1189:    PETSc must be configured with ZLIB package.

1191:    See also examples src/mat/tutorials/ex10.c and src/ksp/ksp/tutorials/ex27.c

1193:    Current HDF5 (MAT-File) limitations:
1194:    This reader currently supports only real MATSEQAIJ, MATMPIAIJ, MATSEQDENSE and MATMPIDENSE matrices.

1196:    Corresponding MatView() is not yet implemented.

1198:    The loaded matrix is actually a transpose of the original one in MATLAB,
1199:    unless you push PETSC_VIEWER_HDF5_MAT format (see examples above).
1200:    With this format, matrix is automatically transposed by PETSc,
1201:    unless the matrix is marked as SPD or symmetric
1202:    (see MatSetOption(), MAT_SPD, MAT_SYMMETRIC).

1204:    References:
1205: 1. MATLAB(R) Documentation, manual page of save(), https://www.mathworks.com/help/matlab/ref/save.html#btox10b-1-version

1207: .seealso: PetscViewerBinaryOpen(), PetscViewerSetType(), MatView(), VecLoad()

1209:  @*/
1210: PetscErrorCode MatLoad(Mat mat,PetscViewer viewer)
1211: {
1213:   PetscBool      flg;


1219:   if (!((PetscObject)mat)->type_name) {
1220:     MatSetType(mat,MATAIJ);
1221:   }

1223:   flg  = PETSC_FALSE;
1224:   PetscOptionsGetBool(((PetscObject)mat)->options,((PetscObject)mat)->prefix,"-matload_symmetric",&flg,NULL);
1225:   if (flg) {
1226:     MatSetOption(mat,MAT_SYMMETRIC,PETSC_TRUE);
1227:     MatSetOption(mat,MAT_SYMMETRY_ETERNAL,PETSC_TRUE);
1228:   }
1229:   flg  = PETSC_FALSE;
1230:   PetscOptionsGetBool(((PetscObject)mat)->options,((PetscObject)mat)->prefix,"-matload_spd",&flg,NULL);
1231:   if (flg) {
1232:     MatSetOption(mat,MAT_SPD,PETSC_TRUE);
1233:   }

1235:   if (!mat->ops->load) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatLoad is not supported for type %s",((PetscObject)mat)->type_name);
1236:   PetscLogEventBegin(MAT_Load,mat,viewer,0,0);
1237:   (*mat->ops->load)(mat,viewer);
1238:   PetscLogEventEnd(MAT_Load,mat,viewer,0,0);
1239:   return(0);
1240: }

1242: static PetscErrorCode MatDestroy_Redundant(Mat_Redundant **redundant)
1243: {
1245:   Mat_Redundant  *redund = *redundant;
1246:   PetscInt       i;

1249:   if (redund){
1250:     if (redund->matseq) { /* via MatCreateSubMatrices()  */
1251:       ISDestroy(&redund->isrow);
1252:       ISDestroy(&redund->iscol);
1253:       MatDestroySubMatrices(1,&redund->matseq);
1254:     } else {
1255:       PetscFree2(redund->send_rank,redund->recv_rank);
1256:       PetscFree(redund->sbuf_j);
1257:       PetscFree(redund->sbuf_a);
1258:       for (i=0; i<redund->nrecvs; i++) {
1259:         PetscFree(redund->rbuf_j[i]);
1260:         PetscFree(redund->rbuf_a[i]);
1261:       }
1262:       PetscFree4(redund->sbuf_nz,redund->rbuf_nz,redund->rbuf_j,redund->rbuf_a);
1263:     }

1265:     if (redund->subcomm) {
1266:       PetscCommDestroy(&redund->subcomm);
1267:     }
1268:     PetscFree(redund);
1269:   }
1270:   return(0);
1271: }

1273: /*@C
1274:    MatDestroy - Frees space taken by a matrix.

1276:    Collective on Mat

1278:    Input Parameter:
1279: .  A - the matrix

1281:    Level: beginner

1283: @*/
1284: PetscErrorCode MatDestroy(Mat *A)
1285: {

1289:   if (!*A) return(0);
1291:   if (--((PetscObject)(*A))->refct > 0) {*A = NULL; return(0);}

1293:   /* if memory was published with SAWs then destroy it */
1294:   PetscObjectSAWsViewOff((PetscObject)*A);
1295:   if ((*A)->ops->destroy) {
1296:     (*(*A)->ops->destroy)(*A);
1297:   }

1299:   PetscFree((*A)->defaultvectype);
1300:   PetscFree((*A)->bsizes);
1301:   PetscFree((*A)->solvertype);
1302:   MatDestroy_Redundant(&(*A)->redundant);
1303:   MatProductClear(*A);
1304:   MatNullSpaceDestroy(&(*A)->nullsp);
1305:   MatNullSpaceDestroy(&(*A)->transnullsp);
1306:   MatNullSpaceDestroy(&(*A)->nearnullsp);
1307:   MatDestroy(&(*A)->schur);
1308:   PetscLayoutDestroy(&(*A)->rmap);
1309:   PetscLayoutDestroy(&(*A)->cmap);
1310:   PetscHeaderDestroy(A);
1311:   return(0);
1312: }

1314: /*@C
1315:    MatSetValues - Inserts or adds a block of values into a matrix.
1316:    These values may be cached, so MatAssemblyBegin() and MatAssemblyEnd()
1317:    MUST be called after all calls to MatSetValues() have been completed.

1319:    Not Collective

1321:    Input Parameters:
1322: +  mat - the matrix
1323: .  v - a logically two-dimensional array of values
1324: .  m, idxm - the number of rows and their global indices
1325: .  n, idxn - the number of columns and their global indices
1326: -  addv - either ADD_VALUES or INSERT_VALUES, where
1327:    ADD_VALUES adds values to any existing entries, and
1328:    INSERT_VALUES replaces existing entries with new values

1330:    Notes:
1331:    If you create the matrix yourself (that is not with a call to DMCreateMatrix()) then you MUST call MatXXXXSetPreallocation() or
1332:       MatSetUp() before using this routine

1334:    By default the values, v, are row-oriented. See MatSetOption() for other options.

1336:    Calls to MatSetValues() with the INSERT_VALUES and ADD_VALUES
1337:    options cannot be mixed without intervening calls to the assembly
1338:    routines.

1340:    MatSetValues() uses 0-based row and column numbers in Fortran
1341:    as well as in C.

1343:    Negative indices may be passed in idxm and idxn, these rows and columns are
1344:    simply ignored. This allows easily inserting element stiffness matrices
1345:    with homogeneous Dirchlet boundary conditions that you don't want represented
1346:    in the matrix.

1348:    Efficiency Alert:
1349:    The routine MatSetValuesBlocked() may offer much better efficiency
1350:    for users of block sparse formats (MATSEQBAIJ and MATMPIBAIJ).

1352:    Level: beginner

1354:    Developer Notes:
1355:     This is labeled with C so does not automatically generate Fortran stubs and interfaces
1356:                     because it requires multiple Fortran interfaces depending on which arguments are scalar or arrays.

1358: .seealso: MatSetOption(), MatAssemblyBegin(), MatAssemblyEnd(), MatSetValuesBlocked(), MatSetValuesLocal(),
1359:           InsertMode, INSERT_VALUES, ADD_VALUES
1360: @*/
1361: PetscErrorCode MatSetValues(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],const PetscScalar v[],InsertMode addv)
1362: {

1368:   if (!m || !n) return(0); /* no values to insert */
1371:   MatCheckPreallocated(mat,1);

1373:   if (mat->insertmode == NOT_SET_VALUES) {
1374:     mat->insertmode = addv;
1375:   } else if (PetscUnlikely(mat->insertmode != addv)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Cannot mix add values and insert values");
1376:   if (PetscDefined(USE_DEBUG)) {
1377:     PetscInt       i,j;

1379:     if (mat->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1380:     if (!mat->ops->setvalues) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);

1382:     for (i=0; i<m; i++) {
1383:       for (j=0; j<n; j++) {
1384:         if (mat->erroriffailure && PetscIsInfOrNanScalar(v[i*n+j]))
1385: #if defined(PETSC_USE_COMPLEX)
1386:           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]);
1387: #else
1388:           SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_FP,"Inserting %g at matrix entry (%D,%D)",(double)v[i*n+j],idxm[i],idxn[j]);
1389: #endif
1390:       }
1391:     }
1392:     for (i=0; i<m; i++) if (idxm[i] >= mat->rmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Cannot insert in row %D, maximum is %D",idxm[i],mat->rmap->N-1);
1393:     for (i=0; i<n; i++) if (idxn[i] >= mat->cmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Cannot insert in column %D, maximum is %D",idxn[i],mat->cmap->N-1);
1394:   }

1396:   if (mat->assembled) {
1397:     mat->was_assembled = PETSC_TRUE;
1398:     mat->assembled     = PETSC_FALSE;
1399:   }
1400:   PetscLogEventBegin(MAT_SetValues,mat,0,0,0);
1401:   (*mat->ops->setvalues)(mat,m,idxm,n,idxn,v,addv);
1402:   PetscLogEventEnd(MAT_SetValues,mat,0,0,0);
1403:   return(0);
1404: }


1407: /*@
1408:    MatSetValuesRowLocal - Inserts a row (block row for BAIJ matrices) of nonzero
1409:         values into a matrix

1411:    Not Collective

1413:    Input Parameters:
1414: +  mat - the matrix
1415: .  row - the (block) row to set
1416: -  v - a logically two-dimensional array of values

1418:    Notes:
1419:    By the values, v, are column-oriented (for the block version) and sorted

1421:    All the nonzeros in the row must be provided

1423:    The matrix must have previously had its column indices set

1425:    The row must belong to this process

1427:    Level: intermediate

1429: .seealso: MatSetOption(), MatAssemblyBegin(), MatAssemblyEnd(), MatSetValuesBlocked(), MatSetValuesLocal(),
1430:           InsertMode, INSERT_VALUES, ADD_VALUES, MatSetValues(), MatSetValuesRow(), MatSetLocalToGlobalMapping()
1431: @*/
1432: PetscErrorCode MatSetValuesRowLocal(Mat mat,PetscInt row,const PetscScalar v[])
1433: {
1435:   PetscInt       globalrow;

1441:   ISLocalToGlobalMappingApply(mat->rmap->mapping,1,&row,&globalrow);
1442:   MatSetValuesRow(mat,globalrow,v);
1443:   return(0);
1444: }

1446: /*@
1447:    MatSetValuesRow - Inserts a row (block row for BAIJ matrices) of nonzero
1448:         values into a matrix

1450:    Not Collective

1452:    Input Parameters:
1453: +  mat - the matrix
1454: .  row - the (block) row to set
1455: -  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

1457:    Notes:
1458:    The values, v, are column-oriented for the block version.

1460:    All the nonzeros in the row must be provided

1462:    THE MATRIX MUST HAVE PREVIOUSLY HAD ITS COLUMN INDICES SET. IT IS RARE THAT THIS ROUTINE IS USED, usually MatSetValues() is used.

1464:    The row must belong to this process

1466:    Level: advanced

1468: .seealso: MatSetOption(), MatAssemblyBegin(), MatAssemblyEnd(), MatSetValuesBlocked(), MatSetValuesLocal(),
1469:           InsertMode, INSERT_VALUES, ADD_VALUES, MatSetValues()
1470: @*/
1471: PetscErrorCode MatSetValuesRow(Mat mat,PetscInt row,const PetscScalar v[])
1472: {

1478:   MatCheckPreallocated(mat,1);
1480:   if (PetscUnlikely(mat->insertmode == ADD_VALUES)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Cannot mix add and insert values");
1481:   if (PetscUnlikely(mat->factortype)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1482:   mat->insertmode = INSERT_VALUES;

1484:   if (mat->assembled) {
1485:     mat->was_assembled = PETSC_TRUE;
1486:     mat->assembled     = PETSC_FALSE;
1487:   }
1488:   PetscLogEventBegin(MAT_SetValues,mat,0,0,0);
1489:   if (!mat->ops->setvaluesrow) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
1490:   (*mat->ops->setvaluesrow)(mat,row,v);
1491:   PetscLogEventEnd(MAT_SetValues,mat,0,0,0);
1492:   return(0);
1493: }

1495: /*@
1496:    MatSetValuesStencil - Inserts or adds a block of values into a matrix.
1497:      Using structured grid indexing

1499:    Not Collective

1501:    Input Parameters:
1502: +  mat - the matrix
1503: .  m - number of rows being entered
1504: .  idxm - grid coordinates (and component number when dof > 1) for matrix rows being entered
1505: .  n - number of columns being entered
1506: .  idxn - grid coordinates (and component number when dof > 1) for matrix columns being entered
1507: .  v - a logically two-dimensional array of values
1508: -  addv - either ADD_VALUES or INSERT_VALUES, where
1509:    ADD_VALUES adds values to any existing entries, and
1510:    INSERT_VALUES replaces existing entries with new values

1512:    Notes:
1513:    By default the values, v, are row-oriented.  See MatSetOption() for other options.

1515:    Calls to MatSetValuesStencil() with the INSERT_VALUES and ADD_VALUES
1516:    options cannot be mixed without intervening calls to the assembly
1517:    routines.

1519:    The grid coordinates are across the entire grid, not just the local portion

1521:    MatSetValuesStencil() uses 0-based row and column numbers in Fortran
1522:    as well as in C.

1524:    For setting/accessing vector values via array coordinates you can use the DMDAVecGetArray() routine

1526:    In order to use this routine you must either obtain the matrix with DMCreateMatrix()
1527:    or call MatSetLocalToGlobalMapping() and MatSetStencil() first.

1529:    The columns and rows in the stencil passed in MUST be contained within the
1530:    ghost region of the given process as set with DMDACreateXXX() or MatSetStencil(). For example,
1531:    if you create a DMDA with an overlap of one grid level and on a particular process its first
1532:    local nonghost x logical coordinate is 6 (so its first ghost x logical coordinate is 5) the
1533:    first i index you can use in your column and row indices in MatSetStencil() is 5.

1535:    In Fortran idxm and idxn should be declared as
1536: $     MatStencil idxm(4,m),idxn(4,n)
1537:    and the values inserted using
1538: $    idxm(MatStencil_i,1) = i
1539: $    idxm(MatStencil_j,1) = j
1540: $    idxm(MatStencil_k,1) = k
1541: $    idxm(MatStencil_c,1) = c
1542:    etc

1544:    For periodic boundary conditions use negative indices for values to the left (below 0; that are to be
1545:    obtained by wrapping values from right edge). For values to the right of the last entry using that index plus one
1546:    etc to obtain values that obtained by wrapping the values from the left edge. This does not work for anything but the
1547:    DM_BOUNDARY_PERIODIC boundary type.

1549:    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
1550:    a single value per point) you can skip filling those indices.

1552:    Inspired by the structured grid interface to the HYPRE package
1553:    (https://computation.llnl.gov/projects/hypre-scalable-linear-solvers-multigrid-methods)

1555:    Efficiency Alert:
1556:    The routine MatSetValuesBlockedStencil() may offer much better efficiency
1557:    for users of block sparse formats (MATSEQBAIJ and MATMPIBAIJ).

1559:    Level: beginner

1561: .seealso: MatSetOption(), MatAssemblyBegin(), MatAssemblyEnd(), MatSetValuesBlocked(), MatSetValuesLocal()
1562:           MatSetValues(), MatSetValuesBlockedStencil(), MatSetStencil(), DMCreateMatrix(), DMDAVecGetArray(), MatStencil
1563: @*/
1564: PetscErrorCode MatSetValuesStencil(Mat mat,PetscInt m,const MatStencil idxm[],PetscInt n,const MatStencil idxn[],const PetscScalar v[],InsertMode addv)
1565: {
1567:   PetscInt       buf[8192],*bufm=NULL,*bufn=NULL,*jdxm,*jdxn;
1568:   PetscInt       j,i,dim = mat->stencil.dim,*dims = mat->stencil.dims+1,tmp;
1569:   PetscInt       *starts = mat->stencil.starts,*dxm = (PetscInt*)idxm,*dxn = (PetscInt*)idxn,sdim = dim - (1 - (PetscInt)mat->stencil.noc);

1572:   if (!m || !n) return(0); /* no values to insert */

1578:   if ((m+n) <= (PetscInt)(sizeof(buf)/sizeof(PetscInt))) {
1579:     jdxm = buf; jdxn = buf+m;
1580:   } else {
1581:     PetscMalloc2(m,&bufm,n,&bufn);
1582:     jdxm = bufm; jdxn = bufn;
1583:   }
1584:   for (i=0; i<m; i++) {
1585:     for (j=0; j<3-sdim; j++) dxm++;
1586:     tmp = *dxm++ - starts[0];
1587:     for (j=0; j<dim-1; j++) {
1588:       if ((*dxm++ - starts[j+1]) < 0 || tmp < 0) tmp = -1;
1589:       else                                       tmp = tmp*dims[j] + *(dxm-1) - starts[j+1];
1590:     }
1591:     if (mat->stencil.noc) dxm++;
1592:     jdxm[i] = tmp;
1593:   }
1594:   for (i=0; i<n; i++) {
1595:     for (j=0; j<3-sdim; j++) dxn++;
1596:     tmp = *dxn++ - starts[0];
1597:     for (j=0; j<dim-1; j++) {
1598:       if ((*dxn++ - starts[j+1]) < 0 || tmp < 0) tmp = -1;
1599:       else                                       tmp = tmp*dims[j] + *(dxn-1) - starts[j+1];
1600:     }
1601:     if (mat->stencil.noc) dxn++;
1602:     jdxn[i] = tmp;
1603:   }
1604:   MatSetValuesLocal(mat,m,jdxm,n,jdxn,v,addv);
1605:   PetscFree2(bufm,bufn);
1606:   return(0);
1607: }

1609: /*@
1610:    MatSetValuesBlockedStencil - Inserts or adds a block of values into a matrix.
1611:      Using structured grid indexing

1613:    Not Collective

1615:    Input Parameters:
1616: +  mat - the matrix
1617: .  m - number of rows being entered
1618: .  idxm - grid coordinates for matrix rows being entered
1619: .  n - number of columns being entered
1620: .  idxn - grid coordinates for matrix columns being entered
1621: .  v - a logically two-dimensional array of values
1622: -  addv - either ADD_VALUES or INSERT_VALUES, where
1623:    ADD_VALUES adds values to any existing entries, and
1624:    INSERT_VALUES replaces existing entries with new values

1626:    Notes:
1627:    By default the values, v, are row-oriented and unsorted.
1628:    See MatSetOption() for other options.

1630:    Calls to MatSetValuesBlockedStencil() with the INSERT_VALUES and ADD_VALUES
1631:    options cannot be mixed without intervening calls to the assembly
1632:    routines.

1634:    The grid coordinates are across the entire grid, not just the local portion

1636:    MatSetValuesBlockedStencil() uses 0-based row and column numbers in Fortran
1637:    as well as in C.

1639:    For setting/accessing vector values via array coordinates you can use the DMDAVecGetArray() routine

1641:    In order to use this routine you must either obtain the matrix with DMCreateMatrix()
1642:    or call MatSetBlockSize(), MatSetLocalToGlobalMapping() and MatSetStencil() first.

1644:    The columns and rows in the stencil passed in MUST be contained within the
1645:    ghost region of the given process as set with DMDACreateXXX() or MatSetStencil(). For example,
1646:    if you create a DMDA with an overlap of one grid level and on a particular process its first
1647:    local nonghost x logical coordinate is 6 (so its first ghost x logical coordinate is 5) the
1648:    first i index you can use in your column and row indices in MatSetStencil() is 5.

1650:    In Fortran idxm and idxn should be declared as
1651: $     MatStencil idxm(4,m),idxn(4,n)
1652:    and the values inserted using
1653: $    idxm(MatStencil_i,1) = i
1654: $    idxm(MatStencil_j,1) = j
1655: $    idxm(MatStencil_k,1) = k
1656:    etc

1658:    Negative indices may be passed in idxm and idxn, these rows and columns are
1659:    simply ignored. This allows easily inserting element stiffness matrices
1660:    with homogeneous Dirchlet boundary conditions that you don't want represented
1661:    in the matrix.

1663:    Inspired by the structured grid interface to the HYPRE package
1664:    (https://computation.llnl.gov/projects/hypre-scalable-linear-solvers-multigrid-methods)

1666:    Level: beginner

1668: .seealso: MatSetOption(), MatAssemblyBegin(), MatAssemblyEnd(), MatSetValuesBlocked(), MatSetValuesLocal()
1669:           MatSetValues(), MatSetValuesStencil(), MatSetStencil(), DMCreateMatrix(), DMDAVecGetArray(), MatStencil,
1670:           MatSetBlockSize(), MatSetLocalToGlobalMapping()
1671: @*/
1672: PetscErrorCode MatSetValuesBlockedStencil(Mat mat,PetscInt m,const MatStencil idxm[],PetscInt n,const MatStencil idxn[],const PetscScalar v[],InsertMode addv)
1673: {
1675:   PetscInt       buf[8192],*bufm=NULL,*bufn=NULL,*jdxm,*jdxn;
1676:   PetscInt       j,i,dim = mat->stencil.dim,*dims = mat->stencil.dims+1,tmp;
1677:   PetscInt       *starts = mat->stencil.starts,*dxm = (PetscInt*)idxm,*dxn = (PetscInt*)idxn,sdim = dim - (1 - (PetscInt)mat->stencil.noc);

1680:   if (!m || !n) return(0); /* no values to insert */

1687:   if ((m+n) <= (PetscInt)(sizeof(buf)/sizeof(PetscInt))) {
1688:     jdxm = buf; jdxn = buf+m;
1689:   } else {
1690:     PetscMalloc2(m,&bufm,n,&bufn);
1691:     jdxm = bufm; jdxn = bufn;
1692:   }
1693:   for (i=0; i<m; i++) {
1694:     for (j=0; j<3-sdim; j++) dxm++;
1695:     tmp = *dxm++ - starts[0];
1696:     for (j=0; j<sdim-1; j++) {
1697:       if ((*dxm++ - starts[j+1]) < 0 || tmp < 0) tmp = -1;
1698:       else                                       tmp = tmp*dims[j] + *(dxm-1) - starts[j+1];
1699:     }
1700:     dxm++;
1701:     jdxm[i] = tmp;
1702:   }
1703:   for (i=0; i<n; i++) {
1704:     for (j=0; j<3-sdim; j++) dxn++;
1705:     tmp = *dxn++ - starts[0];
1706:     for (j=0; j<sdim-1; j++) {
1707:       if ((*dxn++ - starts[j+1]) < 0 || tmp < 0) tmp = -1;
1708:       else                                       tmp = tmp*dims[j] + *(dxn-1) - starts[j+1];
1709:     }
1710:     dxn++;
1711:     jdxn[i] = tmp;
1712:   }
1713:   MatSetValuesBlockedLocal(mat,m,jdxm,n,jdxn,v,addv);
1714:   PetscFree2(bufm,bufn);
1715:   return(0);
1716: }

1718: /*@
1719:    MatSetStencil - Sets the grid information for setting values into a matrix via
1720:         MatSetValuesStencil()

1722:    Not Collective

1724:    Input Parameters:
1725: +  mat - the matrix
1726: .  dim - dimension of the grid 1, 2, or 3
1727: .  dims - number of grid points in x, y, and z direction, including ghost points on your processor
1728: .  starts - starting point of ghost nodes on your processor in x, y, and z direction
1729: -  dof - number of degrees of freedom per node


1732:    Inspired by the structured grid interface to the HYPRE package
1733:    (www.llnl.gov/CASC/hyper)

1735:    For matrices generated with DMCreateMatrix() this routine is automatically called and so not needed by the
1736:    user.

1738:    Level: beginner

1740: .seealso: MatSetOption(), MatAssemblyBegin(), MatAssemblyEnd(), MatSetValuesBlocked(), MatSetValuesLocal()
1741:           MatSetValues(), MatSetValuesBlockedStencil(), MatSetValuesStencil()
1742: @*/
1743: PetscErrorCode MatSetStencil(Mat mat,PetscInt dim,const PetscInt dims[],const PetscInt starts[],PetscInt dof)
1744: {
1745:   PetscInt i;


1752:   mat->stencil.dim = dim + (dof > 1);
1753:   for (i=0; i<dim; i++) {
1754:     mat->stencil.dims[i]   = dims[dim-i-1];      /* copy the values in backwards */
1755:     mat->stencil.starts[i] = starts[dim-i-1];
1756:   }
1757:   mat->stencil.dims[dim]   = dof;
1758:   mat->stencil.starts[dim] = 0;
1759:   mat->stencil.noc         = (PetscBool)(dof == 1);
1760:   return(0);
1761: }

1763: /*@C
1764:    MatSetValuesBlocked - Inserts or adds a block of values into a matrix.

1766:    Not Collective

1768:    Input Parameters:
1769: +  mat - the matrix
1770: .  v - a logically two-dimensional array of values
1771: .  m, idxm - the number of block rows and their global block indices
1772: .  n, idxn - the number of block columns and their global block indices
1773: -  addv - either ADD_VALUES or INSERT_VALUES, where
1774:    ADD_VALUES adds values to any existing entries, and
1775:    INSERT_VALUES replaces existing entries with new values

1777:    Notes:
1778:    If you create the matrix yourself (that is not with a call to DMCreateMatrix()) then you MUST call
1779:    MatXXXXSetPreallocation() or MatSetUp() before using this routine.

1781:    The m and n count the NUMBER of blocks in the row direction and column direction,
1782:    NOT the total number of rows/columns; for example, if the block size is 2 and
1783:    you are passing in values for rows 2,3,4,5  then m would be 2 (not 4).
1784:    The values in idxm would be 1 2; that is the first index for each block divided by
1785:    the block size.

1787:    Note that you must call MatSetBlockSize() when constructing this matrix (before
1788:    preallocating it).

1790:    By default the values, v, are row-oriented, so the layout of
1791:    v is the same as for MatSetValues(). See MatSetOption() for other options.

1793:    Calls to MatSetValuesBlocked() with the INSERT_VALUES and ADD_VALUES
1794:    options cannot be mixed without intervening calls to the assembly
1795:    routines.

1797:    MatSetValuesBlocked() uses 0-based row and column numbers in Fortran
1798:    as well as in C.

1800:    Negative indices may be passed in idxm and idxn, these rows and columns are
1801:    simply ignored. This allows easily inserting element stiffness matrices
1802:    with homogeneous Dirchlet boundary conditions that you don't want represented
1803:    in the matrix.

1805:    Each time an entry is set within a sparse matrix via MatSetValues(),
1806:    internal searching must be done to determine where to place the
1807:    data in the matrix storage space.  By instead inserting blocks of
1808:    entries via MatSetValuesBlocked(), the overhead of matrix assembly is
1809:    reduced.

1811:    Example:
1812: $   Suppose m=n=2 and block size(bs) = 2 The array is
1813: $
1814: $   1  2  | 3  4
1815: $   5  6  | 7  8
1816: $   - - - | - - -
1817: $   9  10 | 11 12
1818: $   13 14 | 15 16
1819: $
1820: $   v[] should be passed in like
1821: $   v[] = [1,2,3,4,5,6,7,8,9,10,11,12,13,14,15,16]
1822: $
1823: $  If you are not using row oriented storage of v (that is you called MatSetOption(mat,MAT_ROW_ORIENTED,PETSC_FALSE)) then
1824: $   v[] = [1,5,9,13,2,6,10,14,3,7,11,15,4,8,12,16]

1826:    Level: intermediate

1828: .seealso: MatSetBlockSize(), MatSetOption(), MatAssemblyBegin(), MatAssemblyEnd(), MatSetValues(), MatSetValuesBlockedLocal()
1829: @*/
1830: PetscErrorCode MatSetValuesBlocked(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],const PetscScalar v[],InsertMode addv)
1831: {

1837:   if (!m || !n) return(0); /* no values to insert */
1841:   MatCheckPreallocated(mat,1);
1842:   if (mat->insertmode == NOT_SET_VALUES) {
1843:     mat->insertmode = addv;
1844:   } else if (PetscUnlikely(mat->insertmode != addv)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Cannot mix add values and insert values");
1845:   if (PetscDefined(USE_DEBUG)) {
1846:     if (mat->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1847:     if (!mat->ops->setvaluesblocked && !mat->ops->setvalues) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
1848:   }
1849:   if (PetscDefined(USE_DEBUG)) {
1850:     PetscInt rbs,cbs,M,N,i;
1851:     MatGetBlockSizes(mat,&rbs,&cbs);
1852:     MatGetSize(mat,&M,&N);
1853:     for (i=0; i<m; i++) {
1854:       if (idxm[i]*rbs >= M) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row block index %D (index %D) greater than row length %D",i,idxm[i],M);
1855:     }
1856:     for (i=0; i<n; i++) {
1857:       if (idxn[i]*cbs >= N) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column block index %D (index %D) great than column length %D",i,idxn[i],N);
1858:     }
1859:   }
1860:   if (mat->assembled) {
1861:     mat->was_assembled = PETSC_TRUE;
1862:     mat->assembled     = PETSC_FALSE;
1863:   }
1864:   PetscLogEventBegin(MAT_SetValues,mat,0,0,0);
1865:   if (mat->ops->setvaluesblocked) {
1866:     (*mat->ops->setvaluesblocked)(mat,m,idxm,n,idxn,v,addv);
1867:   } else {
1868:     PetscInt buf[8192],*bufr=NULL,*bufc=NULL,*iidxm,*iidxn;
1869:     PetscInt i,j,bs,cbs;
1870:     MatGetBlockSizes(mat,&bs,&cbs);
1871:     if (m*bs+n*cbs <= (PetscInt)(sizeof(buf)/sizeof(PetscInt))) {
1872:       iidxm = buf; iidxn = buf + m*bs;
1873:     } else {
1874:       PetscMalloc2(m*bs,&bufr,n*cbs,&bufc);
1875:       iidxm = bufr; iidxn = bufc;
1876:     }
1877:     for (i=0; i<m; i++) {
1878:       for (j=0; j<bs; j++) {
1879:         iidxm[i*bs+j] = bs*idxm[i] + j;
1880:       }
1881:     }
1882:     for (i=0; i<n; i++) {
1883:       for (j=0; j<cbs; j++) {
1884:         iidxn[i*cbs+j] = cbs*idxn[i] + j;
1885:       }
1886:     }
1887:     MatSetValues(mat,m*bs,iidxm,n*cbs,iidxn,v,addv);
1888:     PetscFree2(bufr,bufc);
1889:   }
1890:   PetscLogEventEnd(MAT_SetValues,mat,0,0,0);
1891:   return(0);
1892: }

1894: /*@C
1895:    MatGetValues - Gets a block of values from a matrix.

1897:    Not Collective; can only return values that are owned by the give process

1899:    Input Parameters:
1900: +  mat - the matrix
1901: .  v - a logically two-dimensional array for storing the values
1902: .  m, idxm - the number of rows and their global indices
1903: -  n, idxn - the number of columns and their global indices

1905:    Notes:
1906:      The user must allocate space (m*n PetscScalars) for the values, v.
1907:      The values, v, are then returned in a row-oriented format,
1908:      analogous to that used by default in MatSetValues().

1910:      MatGetValues() uses 0-based row and column numbers in
1911:      Fortran as well as in C.

1913:      MatGetValues() requires that the matrix has been assembled
1914:      with MatAssemblyBegin()/MatAssemblyEnd().  Thus, calls to
1915:      MatSetValues() and MatGetValues() CANNOT be made in succession
1916:      without intermediate matrix assembly.

1918:      Negative row or column indices will be ignored and those locations in v[] will be
1919:      left unchanged.

1921:      For the standard row-based matrix formats, idxm[] can only contain rows owned by the requesting MPI rank.
1922:      That is, rows with global index greater than or equal to restart and less than rend where restart and rend are obtainable
1923:      from MatGetOwnershipRange(mat,&rstart,&rend).

1925:    Level: advanced

1927: .seealso: MatGetRow(), MatCreateSubMatrices(), MatSetValues(), MatGetOwnershipRange(), MatGetValuesLocal()
1928: @*/
1929: PetscErrorCode MatGetValues(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],PetscScalar v[])
1930: {

1936:   if (!m || !n) return(0);
1940:   if (!mat->assembled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
1941:   if (mat->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1942:   if (!mat->ops->getvalues) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
1943:   MatCheckPreallocated(mat,1);

1945:   PetscLogEventBegin(MAT_GetValues,mat,0,0,0);
1946:   (*mat->ops->getvalues)(mat,m,idxm,n,idxn,v);
1947:   PetscLogEventEnd(MAT_GetValues,mat,0,0,0);
1948:   return(0);
1949: }

1951: /*@C
1952:    MatGetValuesLocal - retrieves values from certain locations in a matrix using the local numbering of the indices
1953:      defined previously by MatSetLocalToGlobalMapping()

1955:    Not Collective

1957:    Input Parameters:
1958: +  mat - the matrix
1959: .  nrow, irow - number of rows and their local indices
1960: -  ncol, icol - number of columns and their local indices

1962:    Output Parameter:
1963: .  y -  a logically two-dimensional array of values

1965:    Notes:
1966:      If you create the matrix yourself (that is not with a call to DMCreateMatrix()) then you MUST call MatSetLocalToGlobalMapping() before using this routine.

1968:      This routine can only return values that are owned by the requesting MPI rank. That is, for standard matrix formats, rows that, in the global numbering,
1969:      are greater than or equal to restart and less than rend where restart and rend are obtainable from MatGetOwnershipRange(mat,&rstart,&rend). One can
1970:      determine if the resulting global row associated with the local row r is owned by the requesting MPI rank by applying the ISLocalToGlobalMapping set
1971:      with MatSetLocalToGlobalMapping().

1973:    Developer Notes:
1974:       This is labelled with C so does not automatically generate Fortran stubs and interfaces
1975:       because it requires multiple Fortran interfaces depending on which arguments are scalar or arrays.

1977:    Level: advanced

1979: .seealso:  MatAssemblyBegin(), MatAssemblyEnd(), MatSetValues(), MatSetLocalToGlobalMapping(),
1980:            MatSetValuesLocal(), MatGetValues()
1981: @*/
1982: PetscErrorCode MatGetValuesLocal(Mat mat,PetscInt nrow,const PetscInt irow[],PetscInt ncol,const PetscInt icol[],PetscScalar y[])
1983: {

1989:   MatCheckPreallocated(mat,1);
1990:   if (!nrow || !ncol) return(0); /* no values to retrieve */
1993:   if (PetscDefined(USE_DEBUG)) {
1994:     if (mat->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1995:     if (!mat->ops->getvalueslocal && !mat->ops->getvalues) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
1996:   }
1997:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
1998:   PetscLogEventBegin(MAT_GetValues,mat,0,0,0);
1999:   if (mat->ops->getvalueslocal) {
2000:     (*mat->ops->getvalueslocal)(mat,nrow,irow,ncol,icol,y);
2001:   } else {
2002:     PetscInt buf[8192],*bufr=NULL,*bufc=NULL,*irowm,*icolm;
2003:     if ((nrow+ncol) <= (PetscInt)(sizeof(buf)/sizeof(PetscInt))) {
2004:       irowm = buf; icolm = buf+nrow;
2005:     } else {
2006:       PetscMalloc2(nrow,&bufr,ncol,&bufc);
2007:       irowm = bufr; icolm = bufc;
2008:     }
2009:     if (!mat->rmap->mapping) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"MatGetValuesLocal() cannot proceed without local-to-global row mapping (See MatSetLocalToGlobalMapping()).");
2010:     if (!mat->cmap->mapping) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"MatGetValuesLocal() cannot proceed without local-to-global column mapping (See MatSetLocalToGlobalMapping()).");
2011:     ISLocalToGlobalMappingApply(mat->rmap->mapping,nrow,irow,irowm);
2012:     ISLocalToGlobalMappingApply(mat->cmap->mapping,ncol,icol,icolm);
2013:     MatGetValues(mat,nrow,irowm,ncol,icolm,y);
2014:     PetscFree2(bufr,bufc);
2015:   }
2016:   PetscLogEventEnd(MAT_GetValues,mat,0,0,0);
2017:   return(0);
2018: }

2020: /*@
2021:   MatSetValuesBatch - Adds (ADD_VALUES) many blocks of values into a matrix at once. The blocks must all be square and
2022:   the same size. Currently, this can only be called once and creates the given matrix.

2024:   Not Collective

2026:   Input Parameters:
2027: + mat - the matrix
2028: . nb - the number of blocks
2029: . bs - the number of rows (and columns) in each block
2030: . rows - a concatenation of the rows for each block
2031: - v - a concatenation of logically two-dimensional arrays of values

2033:   Notes:
2034:   In the future, we will extend this routine to handle rectangular blocks, and to allow multiple calls for a given matrix.

2036:   Level: advanced

2038: .seealso: MatSetOption(), MatAssemblyBegin(), MatAssemblyEnd(), MatSetValuesBlocked(), MatSetValuesLocal(),
2039:           InsertMode, INSERT_VALUES, ADD_VALUES, MatSetValues()
2040: @*/
2041: PetscErrorCode MatSetValuesBatch(Mat mat, PetscInt nb, PetscInt bs, PetscInt rows[], const PetscScalar v[])
2042: {

2050:   if (PetscUnlikelyDebug(mat->factortype)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");

2052:   PetscLogEventBegin(MAT_SetValuesBatch,mat,0,0,0);
2053:   if (mat->ops->setvaluesbatch) {
2054:     (*mat->ops->setvaluesbatch)(mat,nb,bs,rows,v);
2055:   } else {
2056:     PetscInt b;
2057:     for (b = 0; b < nb; ++b) {
2058:       MatSetValues(mat, bs, &rows[b*bs], bs, &rows[b*bs], &v[b*bs*bs], ADD_VALUES);
2059:     }
2060:   }
2061:   PetscLogEventEnd(MAT_SetValuesBatch,mat,0,0,0);
2062:   return(0);
2063: }

2065: /*@
2066:    MatSetLocalToGlobalMapping - Sets a local-to-global numbering for use by
2067:    the routine MatSetValuesLocal() to allow users to insert matrix entries
2068:    using a local (per-processor) numbering.

2070:    Not Collective

2072:    Input Parameters:
2073: +  x - the matrix
2074: .  rmapping - row mapping created with ISLocalToGlobalMappingCreate()   or ISLocalToGlobalMappingCreateIS()
2075: - cmapping - column mapping

2077:    Level: intermediate


2080: .seealso:  MatAssemblyBegin(), MatAssemblyEnd(), MatSetValues(), MatSetValuesLocal(), MatGetValuesLocal()
2081: @*/
2082: PetscErrorCode MatSetLocalToGlobalMapping(Mat x,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping)
2083: {


2092:   if (x->ops->setlocaltoglobalmapping) {
2093:     (*x->ops->setlocaltoglobalmapping)(x,rmapping,cmapping);
2094:   } else {
2095:     PetscLayoutSetISLocalToGlobalMapping(x->rmap,rmapping);
2096:     PetscLayoutSetISLocalToGlobalMapping(x->cmap,cmapping);
2097:   }
2098:   return(0);
2099: }


2102: /*@
2103:    MatGetLocalToGlobalMapping - Gets the local-to-global numbering set by MatSetLocalToGlobalMapping()

2105:    Not Collective

2107:    Input Parameters:
2108: .  A - the matrix

2110:    Output Parameters:
2111: + rmapping - row mapping
2112: - cmapping - column mapping

2114:    Level: advanced


2117: .seealso:  MatSetValuesLocal()
2118: @*/
2119: PetscErrorCode MatGetLocalToGlobalMapping(Mat A,ISLocalToGlobalMapping *rmapping,ISLocalToGlobalMapping *cmapping)
2120: {
2126:   if (rmapping) *rmapping = A->rmap->mapping;
2127:   if (cmapping) *cmapping = A->cmap->mapping;
2128:   return(0);
2129: }

2131: /*@
2132:    MatSetLayouts - Sets the PetscLayout objects for rows and columns of a matrix

2134:    Logically Collective on A

2136:    Input Parameters:
2137: +  A - the matrix
2138: . rmap - row layout
2139: - cmap - column layout

2141:    Level: advanced

2143: .seealso:  MatCreateVecs(), MatGetLocalToGlobalMapping(), MatGetLayouts()
2144: @*/
2145: PetscErrorCode MatSetLayouts(Mat A,PetscLayout rmap,PetscLayout cmap)
2146: {


2152:   PetscLayoutReference(rmap,&A->rmap);
2153:   PetscLayoutReference(cmap,&A->cmap);
2154:   return(0);
2155: }

2157: /*@
2158:    MatGetLayouts - Gets the PetscLayout objects for rows and columns

2160:    Not Collective

2162:    Input Parameters:
2163: .  A - the matrix

2165:    Output Parameters:
2166: + rmap - row layout
2167: - cmap - column layout

2169:    Level: advanced

2171: .seealso:  MatCreateVecs(), MatGetLocalToGlobalMapping(), MatSetLayouts()
2172: @*/
2173: PetscErrorCode MatGetLayouts(Mat A,PetscLayout *rmap,PetscLayout *cmap)
2174: {
2180:   if (rmap) *rmap = A->rmap;
2181:   if (cmap) *cmap = A->cmap;
2182:   return(0);
2183: }

2185: /*@C
2186:    MatSetValuesLocal - Inserts or adds values into certain locations of a matrix,
2187:    using a local numbering of the nodes.

2189:    Not Collective

2191:    Input Parameters:
2192: +  mat - the matrix
2193: .  nrow, irow - number of rows and their local indices
2194: .  ncol, icol - number of columns and their local indices
2195: .  y -  a logically two-dimensional array of values
2196: -  addv - either INSERT_VALUES or ADD_VALUES, where
2197:    ADD_VALUES adds values to any existing entries, and
2198:    INSERT_VALUES replaces existing entries with new values

2200:    Notes:
2201:    If you create the matrix yourself (that is not with a call to DMCreateMatrix()) then you MUST call MatXXXXSetPreallocation() or
2202:       MatSetUp() before using this routine

2204:    If you create the matrix yourself (that is not with a call to DMCreateMatrix()) then you MUST call MatSetLocalToGlobalMapping() before using this routine

2206:    Calls to MatSetValuesLocal() with the INSERT_VALUES and ADD_VALUES
2207:    options cannot be mixed without intervening calls to the assembly
2208:    routines.

2210:    These values may be cached, so MatAssemblyBegin() and MatAssemblyEnd()
2211:    MUST be called after all calls to MatSetValuesLocal() have been completed.

2213:    Level: intermediate

2215:    Developer Notes:
2216:     This is labeled with C so does not automatically generate Fortran stubs and interfaces
2217:                     because it requires multiple Fortran interfaces depending on which arguments are scalar or arrays.

2219: .seealso:  MatAssemblyBegin(), MatAssemblyEnd(), MatSetValues(), MatSetLocalToGlobalMapping(),
2220:            MatSetValueLocal(), MatGetValuesLocal()
2221: @*/
2222: PetscErrorCode MatSetValuesLocal(Mat mat,PetscInt nrow,const PetscInt irow[],PetscInt ncol,const PetscInt icol[],const PetscScalar y[],InsertMode addv)
2223: {

2229:   MatCheckPreallocated(mat,1);
2230:   if (!nrow || !ncol) return(0); /* no values to insert */
2233:   if (mat->insertmode == NOT_SET_VALUES) {
2234:     mat->insertmode = addv;
2235:   }
2236:   else if (PetscUnlikely(mat->insertmode != addv)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Cannot mix add values and insert values");
2237:   if (PetscDefined(USE_DEBUG)) {
2238:     if (mat->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2239:     if (!mat->ops->setvalueslocal && !mat->ops->setvalues) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
2240:   }

2242:   if (mat->assembled) {
2243:     mat->was_assembled = PETSC_TRUE;
2244:     mat->assembled     = PETSC_FALSE;
2245:   }
2246:   PetscLogEventBegin(MAT_SetValues,mat,0,0,0);
2247:   if (mat->ops->setvalueslocal) {
2248:     (*mat->ops->setvalueslocal)(mat,nrow,irow,ncol,icol,y,addv);
2249:   } else {
2250:     PetscInt buf[8192],*bufr=NULL,*bufc=NULL,*irowm,*icolm;
2251:     if ((nrow+ncol) <= (PetscInt)(sizeof(buf)/sizeof(PetscInt))) {
2252:       irowm = buf; icolm = buf+nrow;
2253:     } else {
2254:       PetscMalloc2(nrow,&bufr,ncol,&bufc);
2255:       irowm = bufr; icolm = bufc;
2256:     }
2257:     if (!mat->rmap->mapping) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"MatSetValuesLocal() cannot proceed without local-to-global row mapping (See MatSetLocalToGlobalMapping()).");
2258:     if (!mat->cmap->mapping) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"MatSetValuesLocal() cannot proceed without local-to-global column mapping (See MatSetLocalToGlobalMapping()).");
2259:     ISLocalToGlobalMappingApply(mat->rmap->mapping,nrow,irow,irowm);
2260:     ISLocalToGlobalMappingApply(mat->cmap->mapping,ncol,icol,icolm);
2261:     MatSetValues(mat,nrow,irowm,ncol,icolm,y,addv);
2262:     PetscFree2(bufr,bufc);
2263:   }
2264:   PetscLogEventEnd(MAT_SetValues,mat,0,0,0);
2265:   return(0);
2266: }

2268: /*@C
2269:    MatSetValuesBlockedLocal - Inserts or adds values into certain locations of a matrix,
2270:    using a local ordering of the nodes a block at a time.

2272:    Not Collective

2274:    Input Parameters:
2275: +  x - the matrix
2276: .  nrow, irow - number of rows and their local indices
2277: .  ncol, icol - number of columns and their local indices
2278: .  y -  a logically two-dimensional array of values
2279: -  addv - either INSERT_VALUES or ADD_VALUES, where
2280:    ADD_VALUES adds values to any existing entries, and
2281:    INSERT_VALUES replaces existing entries with new values

2283:    Notes:
2284:    If you create the matrix yourself (that is not with a call to DMCreateMatrix()) then you MUST call MatXXXXSetPreallocation() or
2285:       MatSetUp() before using this routine

2287:    If you create the matrix yourself (that is not with a call to DMCreateMatrix()) then you MUST call MatSetBlockSize() and MatSetLocalToGlobalMapping()
2288:       before using this routineBefore calling MatSetValuesLocal(), the user must first set the

2290:    Calls to MatSetValuesBlockedLocal() with the INSERT_VALUES and ADD_VALUES
2291:    options cannot be mixed without intervening calls to the assembly
2292:    routines.

2294:    These values may be cached, so MatAssemblyBegin() and MatAssemblyEnd()
2295:    MUST be called after all calls to MatSetValuesBlockedLocal() have been completed.

2297:    Level: intermediate

2299:    Developer Notes:
2300:     This is labeled with C so does not automatically generate Fortran stubs and interfaces
2301:                     because it requires multiple Fortran interfaces depending on which arguments are scalar or arrays.

2303: .seealso:  MatSetBlockSize(), MatSetLocalToGlobalMapping(), MatAssemblyBegin(), MatAssemblyEnd(),
2304:            MatSetValuesLocal(),  MatSetValuesBlocked()
2305: @*/
2306: PetscErrorCode MatSetValuesBlockedLocal(Mat mat,PetscInt nrow,const PetscInt irow[],PetscInt ncol,const PetscInt icol[],const PetscScalar y[],InsertMode addv)
2307: {

2313:   MatCheckPreallocated(mat,1);
2314:   if (!nrow || !ncol) return(0); /* no values to insert */
2318:   if (mat->insertmode == NOT_SET_VALUES) {
2319:     mat->insertmode = addv;
2320:   } else if (PetscUnlikely(mat->insertmode != addv)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Cannot mix add values and insert values");
2321:   if (PetscDefined(USE_DEBUG)) {
2322:     if (mat->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2323:     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);
2324:   }

2326:   if (mat->assembled) {
2327:     mat->was_assembled = PETSC_TRUE;
2328:     mat->assembled     = PETSC_FALSE;
2329:   }
2330:   if (PetscUnlikelyDebug(mat->rmap->mapping)) { /* Condition on the mapping existing, because MatSetValuesBlockedLocal_IS does not require it to be set. */
2331:     PetscInt irbs, rbs;
2332:     MatGetBlockSizes(mat, &rbs, NULL);
2333:     ISLocalToGlobalMappingGetBlockSize(mat->rmap->mapping,&irbs);
2334:     if (rbs != irbs) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Different row block sizes! mat %D, row l2g map %D",rbs,irbs);
2335:   }
2336:   if (PetscUnlikelyDebug(mat->cmap->mapping)) {
2337:     PetscInt icbs, cbs;
2338:     MatGetBlockSizes(mat,NULL,&cbs);
2339:     ISLocalToGlobalMappingGetBlockSize(mat->cmap->mapping,&icbs);
2340:     if (cbs != icbs) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Different col block sizes! mat %D, col l2g map %D",cbs,icbs);
2341:   }
2342:   PetscLogEventBegin(MAT_SetValues,mat,0,0,0);
2343:   if (mat->ops->setvaluesblockedlocal) {
2344:     (*mat->ops->setvaluesblockedlocal)(mat,nrow,irow,ncol,icol,y,addv);
2345:   } else {
2346:     PetscInt buf[8192],*bufr=NULL,*bufc=NULL,*irowm,*icolm;
2347:     if ((nrow+ncol) <= (PetscInt)(sizeof(buf)/sizeof(PetscInt))) {
2348:       irowm = buf; icolm = buf + nrow;
2349:     } else {
2350:       PetscMalloc2(nrow,&bufr,ncol,&bufc);
2351:       irowm = bufr; icolm = bufc;
2352:     }
2353:     ISLocalToGlobalMappingApplyBlock(mat->rmap->mapping,nrow,irow,irowm);
2354:     ISLocalToGlobalMappingApplyBlock(mat->cmap->mapping,ncol,icol,icolm);
2355:     MatSetValuesBlocked(mat,nrow,irowm,ncol,icolm,y,addv);
2356:     PetscFree2(bufr,bufc);
2357:   }
2358:   PetscLogEventEnd(MAT_SetValues,mat,0,0,0);
2359:   return(0);
2360: }

2362: /*@
2363:    MatMultDiagonalBlock - Computes the matrix-vector product, y = Dx. Where D is defined by the inode or block structure of the diagonal

2365:    Collective on Mat

2367:    Input Parameters:
2368: +  mat - the matrix
2369: -  x   - the vector to be multiplied

2371:    Output Parameters:
2372: .  y - the result

2374:    Notes:
2375:    The vectors x and y cannot be the same.  I.e., one cannot
2376:    call MatMult(A,y,y).

2378:    Level: developer

2380: .seealso: MatMultTranspose(), MatMultAdd(), MatMultTransposeAdd()
2381: @*/
2382: PetscErrorCode MatMultDiagonalBlock(Mat mat,Vec x,Vec y)
2383: {


2392:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2393:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2394:   if (x == y) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"x and y must be different vectors");
2395:   MatCheckPreallocated(mat,1);

2397:   if (!mat->ops->multdiagonalblock) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Matrix type %s does not have a multiply defined",((PetscObject)mat)->type_name);
2398:   (*mat->ops->multdiagonalblock)(mat,x,y);
2399:   PetscObjectStateIncrease((PetscObject)y);
2400:   return(0);
2401: }

2403: /* --------------------------------------------------------*/
2404: /*@
2405:    MatMult - Computes the matrix-vector product, y = Ax.

2407:    Neighbor-wise Collective on Mat

2409:    Input Parameters:
2410: +  mat - the matrix
2411: -  x   - the vector to be multiplied

2413:    Output Parameters:
2414: .  y - the result

2416:    Notes:
2417:    The vectors x and y cannot be the same.  I.e., one cannot
2418:    call MatMult(A,y,y).

2420:    Level: beginner

2422: .seealso: MatMultTranspose(), MatMultAdd(), MatMultTransposeAdd()
2423: @*/
2424: PetscErrorCode MatMult(Mat mat,Vec x,Vec y)
2425: {

2433:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2434:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2435:   if (x == y) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"x and y must be different vectors");
2436:   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);
2437:   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);
2438:   if (mat->cmap->n != x->map->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat mat,Vec x: local dim %D %D",mat->cmap->n,x->map->n);
2439:   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);
2440:   VecSetErrorIfLocked(y,3);
2441:   if (mat->erroriffailure) {VecValidValues(x,2,PETSC_TRUE);}
2442:   MatCheckPreallocated(mat,1);

2444:   VecLockReadPush(x);
2445:   if (!mat->ops->mult) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Matrix type %s does not have a multiply defined",((PetscObject)mat)->type_name);
2446:   PetscLogEventBegin(MAT_Mult,mat,x,y,0);
2447:   (*mat->ops->mult)(mat,x,y);
2448:   PetscLogEventEnd(MAT_Mult,mat,x,y,0);
2449:   if (mat->erroriffailure) {VecValidValues(y,3,PETSC_FALSE);}
2450:   VecLockReadPop(x);
2451:   return(0);
2452: }

2454: /*@
2455:    MatMultTranspose - Computes matrix transpose times a vector y = A^T * x.

2457:    Neighbor-wise Collective on Mat

2459:    Input Parameters:
2460: +  mat - the matrix
2461: -  x   - the vector to be multiplied

2463:    Output Parameters:
2464: .  y - the result

2466:    Notes:
2467:    The vectors x and y cannot be the same.  I.e., one cannot
2468:    call MatMultTranspose(A,y,y).

2470:    For complex numbers this does NOT compute the Hermitian (complex conjugate) transpose multiple,
2471:    use MatMultHermitianTranspose()

2473:    Level: beginner

2475: .seealso: MatMult(), MatMultAdd(), MatMultTransposeAdd(), MatMultHermitianTranspose(), MatTranspose()
2476: @*/
2477: PetscErrorCode MatMultTranspose(Mat mat,Vec x,Vec y)
2478: {
2479:   PetscErrorCode (*op)(Mat,Vec,Vec)=NULL,ierr;


2487:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2488:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2489:   if (x == y) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"x and y must be different vectors");
2490:   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);
2491:   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);
2492:   if (mat->cmap->n != y->map->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat mat,Vec y: local dim %D %D",mat->cmap->n,y->map->n);
2493:   if (mat->rmap->n != x->map->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat mat,Vec x: local dim %D %D",mat->rmap->n,x->map->n);
2494:   if (mat->erroriffailure) {VecValidValues(x,2,PETSC_TRUE);}
2495:   MatCheckPreallocated(mat,1);

2497:   if (!mat->ops->multtranspose) {
2498:     if (mat->symmetric && mat->ops->mult) op = mat->ops->mult;
2499:     if (!op) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Matrix type %s does not have a multiply transpose defined or is symmetric and does not have a multiply defined",((PetscObject)mat)->type_name);
2500:   } else op = mat->ops->multtranspose;
2501:   PetscLogEventBegin(MAT_MultTranspose,mat,x,y,0);
2502:   VecLockReadPush(x);
2503:   (*op)(mat,x,y);
2504:   VecLockReadPop(x);
2505:   PetscLogEventEnd(MAT_MultTranspose,mat,x,y,0);
2506:   PetscObjectStateIncrease((PetscObject)y);
2507:   if (mat->erroriffailure) {VecValidValues(y,3,PETSC_FALSE);}
2508:   return(0);
2509: }

2511: /*@
2512:    MatMultHermitianTranspose - Computes matrix Hermitian transpose times a vector.

2514:    Neighbor-wise Collective on Mat

2516:    Input Parameters:
2517: +  mat - the matrix
2518: -  x   - the vector to be multilplied

2520:    Output Parameters:
2521: .  y - the result

2523:    Notes:
2524:    The vectors x and y cannot be the same.  I.e., one cannot
2525:    call MatMultHermitianTranspose(A,y,y).

2527:    Also called the conjugate transpose, complex conjugate transpose, or adjoint.

2529:    For real numbers MatMultTranspose() and MatMultHermitianTranspose() are identical.

2531:    Level: beginner

2533: .seealso: MatMult(), MatMultAdd(), MatMultHermitianTransposeAdd(), MatMultTranspose()
2534: @*/
2535: PetscErrorCode MatMultHermitianTranspose(Mat mat,Vec x,Vec y)
2536: {


2545:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2546:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2547:   if (x == y) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"x and y must be different vectors");
2548:   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);
2549:   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);
2550:   if (mat->cmap->n != y->map->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat mat,Vec y: local dim %D %D",mat->cmap->n,y->map->n);
2551:   if (mat->rmap->n != x->map->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat mat,Vec x: local dim %D %D",mat->rmap->n,x->map->n);
2552:   MatCheckPreallocated(mat,1);

2554:   PetscLogEventBegin(MAT_MultHermitianTranspose,mat,x,y,0);
2555: #if defined(PETSC_USE_COMPLEX)
2556:   if (mat->ops->multhermitiantranspose || (mat->hermitian && mat->ops->mult)) {
2557:     VecLockReadPush(x);
2558:     if (mat->ops->multhermitiantranspose) {
2559:       (*mat->ops->multhermitiantranspose)(mat,x,y);
2560:     } else {
2561:       (*mat->ops->mult)(mat,x,y);
2562:     }
2563:     VecLockReadPop(x);
2564:   } else {
2565:     Vec w;
2566:     VecDuplicate(x,&w);
2567:     VecCopy(x,w);
2568:     VecConjugate(w);
2569:     MatMultTranspose(mat,w,y);
2570:     VecDestroy(&w);
2571:     VecConjugate(y);
2572:   }
2573:   PetscObjectStateIncrease((PetscObject)y);
2574: #else
2575:   MatMultTranspose(mat,x,y);
2576: #endif
2577:   PetscLogEventEnd(MAT_MultHermitianTranspose,mat,x,y,0);
2578:   return(0);
2579: }

2581: /*@
2582:     MatMultAdd -  Computes v3 = v2 + A * v1.

2584:     Neighbor-wise Collective on Mat

2586:     Input Parameters:
2587: +   mat - the matrix
2588: -   v1, v2 - the vectors

2590:     Output Parameters:
2591: .   v3 - the result

2593:     Notes:
2594:     The vectors v1 and v3 cannot be the same.  I.e., one cannot
2595:     call MatMultAdd(A,v1,v2,v1).

2597:     Level: beginner

2599: .seealso: MatMultTranspose(), MatMult(), MatMultTransposeAdd()
2600: @*/
2601: PetscErrorCode MatMultAdd(Mat mat,Vec v1,Vec v2,Vec v3)
2602: {


2612:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2613:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2614:   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);
2615:   /* 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);
2616:      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); */
2617:   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);
2618:   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);
2619:   if (v1 == v3) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_IDN,"v1 and v3 must be different vectors");
2620:   MatCheckPreallocated(mat,1);

2622:   if (!mat->ops->multadd) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"No MatMultAdd() for matrix type %s",((PetscObject)mat)->type_name);
2623:   PetscLogEventBegin(MAT_MultAdd,mat,v1,v2,v3);
2624:   VecLockReadPush(v1);
2625:   (*mat->ops->multadd)(mat,v1,v2,v3);
2626:   VecLockReadPop(v1);
2627:   PetscLogEventEnd(MAT_MultAdd,mat,v1,v2,v3);
2628:   PetscObjectStateIncrease((PetscObject)v3);
2629:   return(0);
2630: }

2632: /*@
2633:    MatMultTransposeAdd - Computes v3 = v2 + A' * v1.

2635:    Neighbor-wise Collective on Mat

2637:    Input Parameters:
2638: +  mat - the matrix
2639: -  v1, v2 - the vectors

2641:    Output Parameters:
2642: .  v3 - the result

2644:    Notes:
2645:    The vectors v1 and v3 cannot be the same.  I.e., one cannot
2646:    call MatMultTransposeAdd(A,v1,v2,v1).

2648:    Level: beginner

2650: .seealso: MatMultTranspose(), MatMultAdd(), MatMult()
2651: @*/
2652: PetscErrorCode MatMultTransposeAdd(Mat mat,Vec v1,Vec v2,Vec v3)
2653: {


2663:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2664:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2665:   if (!mat->ops->multtransposeadd) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
2666:   if (v1 == v3) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_IDN,"v1 and v3 must be different vectors");
2667:   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);
2668:   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);
2669:   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);
2670:   MatCheckPreallocated(mat,1);

2672:   PetscLogEventBegin(MAT_MultTransposeAdd,mat,v1,v2,v3);
2673:   VecLockReadPush(v1);
2674:   (*mat->ops->multtransposeadd)(mat,v1,v2,v3);
2675:   VecLockReadPop(v1);
2676:   PetscLogEventEnd(MAT_MultTransposeAdd,mat,v1,v2,v3);
2677:   PetscObjectStateIncrease((PetscObject)v3);
2678:   return(0);
2679: }

2681: /*@
2682:    MatMultHermitianTransposeAdd - Computes v3 = v2 + A^H * v1.

2684:    Neighbor-wise Collective on Mat

2686:    Input Parameters:
2687: +  mat - the matrix
2688: -  v1, v2 - the vectors

2690:    Output Parameters:
2691: .  v3 - the result

2693:    Notes:
2694:    The vectors v1 and v3 cannot be the same.  I.e., one cannot
2695:    call MatMultHermitianTransposeAdd(A,v1,v2,v1).

2697:    Level: beginner

2699: .seealso: MatMultHermitianTranspose(), MatMultTranspose(), MatMultAdd(), MatMult()
2700: @*/
2701: PetscErrorCode MatMultHermitianTransposeAdd(Mat mat,Vec v1,Vec v2,Vec v3)
2702: {


2712:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2713:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2714:   if (v1 == v3) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_IDN,"v1 and v3 must be different vectors");
2715:   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);
2716:   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);
2717:   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);
2718:   MatCheckPreallocated(mat,1);

2720:   PetscLogEventBegin(MAT_MultHermitianTransposeAdd,mat,v1,v2,v3);
2721:   VecLockReadPush(v1);
2722:   if (mat->ops->multhermitiantransposeadd) {
2723:     (*mat->ops->multhermitiantransposeadd)(mat,v1,v2,v3);
2724:   } else {
2725:     Vec w,z;
2726:     VecDuplicate(v1,&w);
2727:     VecCopy(v1,w);
2728:     VecConjugate(w);
2729:     VecDuplicate(v3,&z);
2730:     MatMultTranspose(mat,w,z);
2731:     VecDestroy(&w);
2732:     VecConjugate(z);
2733:     if (v2 != v3) {
2734:       VecWAXPY(v3,1.0,v2,z);
2735:     } else {
2736:       VecAXPY(v3,1.0,z);
2737:     }
2738:     VecDestroy(&z);
2739:   }
2740:   VecLockReadPop(v1);
2741:   PetscLogEventEnd(MAT_MultHermitianTransposeAdd,mat,v1,v2,v3);
2742:   PetscObjectStateIncrease((PetscObject)v3);
2743:   return(0);
2744: }

2746: /*@
2747:    MatMultConstrained - The inner multiplication routine for a
2748:    constrained matrix P^T A P.

2750:    Neighbor-wise Collective on Mat

2752:    Input Parameters:
2753: +  mat - the matrix
2754: -  x   - the vector to be multilplied

2756:    Output Parameters:
2757: .  y - the result

2759:    Notes:
2760:    The vectors x and y cannot be the same.  I.e., one cannot
2761:    call MatMult(A,y,y).

2763:    Level: beginner

2765: .seealso: MatMult(), MatMultTranspose(), MatMultAdd(), MatMultTransposeAdd()
2766: @*/
2767: PetscErrorCode MatMultConstrained(Mat mat,Vec x,Vec y)
2768: {

2775:   if (!mat->assembled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2776:   if (mat->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2777:   if (x == y) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"x and y must be different vectors");
2778:   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);
2779:   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);
2780:   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);

2782:   PetscLogEventBegin(MAT_MultConstrained,mat,x,y,0);
2783:   VecLockReadPush(x);
2784:   (*mat->ops->multconstrained)(mat,x,y);
2785:   VecLockReadPop(x);
2786:   PetscLogEventEnd(MAT_MultConstrained,mat,x,y,0);
2787:   PetscObjectStateIncrease((PetscObject)y);
2788:   return(0);
2789: }

2791: /*@
2792:    MatMultTransposeConstrained - The inner multiplication routine for a
2793:    constrained matrix P^T A^T P.

2795:    Neighbor-wise Collective on Mat

2797:    Input Parameters:
2798: +  mat - the matrix
2799: -  x   - the vector to be multilplied

2801:    Output Parameters:
2802: .  y - the result

2804:    Notes:
2805:    The vectors x and y cannot be the same.  I.e., one cannot
2806:    call MatMult(A,y,y).

2808:    Level: beginner

2810: .seealso: MatMult(), MatMultTranspose(), MatMultAdd(), MatMultTransposeAdd()
2811: @*/
2812: PetscErrorCode MatMultTransposeConstrained(Mat mat,Vec x,Vec y)
2813: {

2820:   if (!mat->assembled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2821:   if (mat->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2822:   if (x == y) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"x and y must be different vectors");
2823:   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);
2824:   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);

2826:   PetscLogEventBegin(MAT_MultConstrained,mat,x,y,0);
2827:   (*mat->ops->multtransposeconstrained)(mat,x,y);
2828:   PetscLogEventEnd(MAT_MultConstrained,mat,x,y,0);
2829:   PetscObjectStateIncrease((PetscObject)y);
2830:   return(0);
2831: }

2833: /*@C
2834:    MatGetFactorType - gets the type of factorization it is

2836:    Not Collective

2838:    Input Parameters:
2839: .  mat - the matrix

2841:    Output Parameters:
2842: .  t - the type, one of MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC,MAT_FACTOR_ILUDT

2844:    Level: intermediate

2846: .seealso: MatFactorType, MatGetFactor(), MatSetFactorType()
2847: @*/
2848: PetscErrorCode MatGetFactorType(Mat mat,MatFactorType *t)
2849: {
2854:   *t = mat->factortype;
2855:   return(0);
2856: }

2858: /*@C
2859:    MatSetFactorType - sets the type of factorization it is

2861:    Logically Collective on Mat

2863:    Input Parameters:
2864: +  mat - the matrix
2865: -  t - the type, one of MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC,MAT_FACTOR_ILUDT

2867:    Level: intermediate

2869: .seealso: MatFactorType, MatGetFactor(), MatGetFactorType()
2870: @*/
2871: PetscErrorCode MatSetFactorType(Mat mat, MatFactorType t)
2872: {
2876:   mat->factortype = t;
2877:   return(0);
2878: }

2880: /* ------------------------------------------------------------*/
2881: /*@C
2882:    MatGetInfo - Returns information about matrix storage (number of
2883:    nonzeros, memory, etc.).

2885:    Collective on Mat if MAT_GLOBAL_MAX or MAT_GLOBAL_SUM is used as the flag

2887:    Input Parameters:
2888: .  mat - the matrix

2890:    Output Parameters:
2891: +  flag - flag indicating the type of parameters to be returned
2892:    (MAT_LOCAL - local matrix, MAT_GLOBAL_MAX - maximum over all processors,
2893:    MAT_GLOBAL_SUM - sum over all processors)
2894: -  info - matrix information context

2896:    Notes:
2897:    The MatInfo context contains a variety of matrix data, including
2898:    number of nonzeros allocated and used, number of mallocs during
2899:    matrix assembly, etc.  Additional information for factored matrices
2900:    is provided (such as the fill ratio, number of mallocs during
2901:    factorization, etc.).  Much of this info is printed to PETSC_STDOUT
2902:    when using the runtime options
2903: $       -info -mat_view ::ascii_info

2905:    Example for C/C++ Users:
2906:    See the file ${PETSC_DIR}/include/petscmat.h for a complete list of
2907:    data within the MatInfo context.  For example,
2908: .vb
2909:       MatInfo info;
2910:       Mat     A;
2911:       double  mal, nz_a, nz_u;

2913:       MatGetInfo(A,MAT_LOCAL,&info);
2914:       mal  = info.mallocs;
2915:       nz_a = info.nz_allocated;
2916: .ve

2918:    Example for Fortran Users:
2919:    Fortran users should declare info as a double precision
2920:    array of dimension MAT_INFO_SIZE, and then extract the parameters
2921:    of interest.  See the file ${PETSC_DIR}/include/petsc/finclude/petscmat.h
2922:    a complete list of parameter names.
2923: .vb
2924:       double  precision info(MAT_INFO_SIZE)
2925:       double  precision mal, nz_a
2926:       Mat     A
2927:       integer ierr

2929:       call MatGetInfo(A,MAT_LOCAL,info,ierr)
2930:       mal = info(MAT_INFO_MALLOCS)
2931:       nz_a = info(MAT_INFO_NZ_ALLOCATED)
2932: .ve

2934:     Level: intermediate

2936:     Developer Note: fortran interface is not autogenerated as the f90
2937:     interface defintion cannot be generated correctly [due to MatInfo]

2939: .seealso: MatStashGetInfo()

2941: @*/
2942: PetscErrorCode MatGetInfo(Mat mat,MatInfoType flag,MatInfo *info)
2943: {

2950:   if (!mat->ops->getinfo) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
2951:   MatCheckPreallocated(mat,1);
2952:   (*mat->ops->getinfo)(mat,flag,info);
2953:   return(0);
2954: }

2956: /*
2957:    This is used by external packages where it is not easy to get the info from the actual
2958:    matrix factorization.
2959: */
2960: PetscErrorCode MatGetInfo_External(Mat A,MatInfoType flag,MatInfo *info)
2961: {

2965:   PetscMemzero(info,sizeof(MatInfo));
2966:   return(0);
2967: }

2969: /* ----------------------------------------------------------*/

2971: /*@C
2972:    MatLUFactor - Performs in-place LU factorization of matrix.

2974:    Collective on Mat

2976:    Input Parameters:
2977: +  mat - the matrix
2978: .  row - row permutation
2979: .  col - column permutation
2980: -  info - options for factorization, includes
2981: $          fill - expected fill as ratio of original fill.
2982: $          dtcol - pivot tolerance (0 no pivot, 1 full column pivoting)
2983: $                   Run with the option -info to determine an optimal value to use

2985:    Notes:
2986:    Most users should employ the simplified KSP interface for linear solvers
2987:    instead of working directly with matrix algebra routines such as this.
2988:    See, e.g., KSPCreate().

2990:    This changes the state of the matrix to a factored matrix; it cannot be used
2991:    for example with MatSetValues() unless one first calls MatSetUnfactored().

2993:    Level: developer

2995: .seealso: MatLUFactorSymbolic(), MatLUFactorNumeric(), MatCholeskyFactor(),
2996:           MatGetOrdering(), MatSetUnfactored(), MatFactorInfo, MatGetFactor()

2998:     Developer Note: fortran interface is not autogenerated as the f90
2999:     interface defintion cannot be generated correctly [due to MatFactorInfo]

3001: @*/
3002: PetscErrorCode MatLUFactor(Mat mat,IS row,IS col,const MatFactorInfo *info)
3003: {
3005:   MatFactorInfo  tinfo;

3013:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
3014:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
3015:   if (!mat->ops->lufactor) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
3016:   MatCheckPreallocated(mat,1);
3017:   if (!info) {
3018:     MatFactorInfoInitialize(&tinfo);
3019:     info = &tinfo;
3020:   }

3022:   PetscLogEventBegin(MAT_LUFactor,mat,row,col,0);
3023:   (*mat->ops->lufactor)(mat,row,col,info);
3024:   PetscLogEventEnd(MAT_LUFactor,mat,row,col,0);
3025:   PetscObjectStateIncrease((PetscObject)mat);
3026:   return(0);
3027: }

3029: /*@C
3030:    MatILUFactor - Performs in-place ILU factorization of matrix.

3032:    Collective on Mat

3034:    Input Parameters:
3035: +  mat - the matrix
3036: .  row - row permutation
3037: .  col - column permutation
3038: -  info - structure containing
3039: $      levels - number of levels of fill.
3040: $      expected fill - as ratio of original fill.
3041: $      1 or 0 - indicating force fill on diagonal (improves robustness for matrices
3042:                 missing diagonal entries)

3044:    Notes:
3045:    Probably really in-place only when level of fill is zero, otherwise allocates
3046:    new space to store factored matrix and deletes previous memory.

3048:    Most users should employ the simplified KSP interface for linear solvers
3049:    instead of working directly with matrix algebra routines such as this.
3050:    See, e.g., KSPCreate().

3052:    Level: developer

3054: .seealso: MatILUFactorSymbolic(), MatLUFactorNumeric(), MatCholeskyFactor(), MatFactorInfo

3056:     Developer Note: fortran interface is not autogenerated as the f90
3057:     interface defintion cannot be generated correctly [due to MatFactorInfo]

3059: @*/
3060: PetscErrorCode MatILUFactor(Mat mat,IS row,IS col,const MatFactorInfo *info)
3061: {

3070:   if (mat->rmap->N != mat->cmap->N) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONG,"matrix must be square");
3071:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
3072:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
3073:   if (!mat->ops->ilufactor) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
3074:   MatCheckPreallocated(mat,1);

3076:   PetscLogEventBegin(MAT_ILUFactor,mat,row,col,0);
3077:   (*mat->ops->ilufactor)(mat,row,col,info);
3078:   PetscLogEventEnd(MAT_ILUFactor,mat,row,col,0);
3079:   PetscObjectStateIncrease((PetscObject)mat);
3080:   return(0);
3081: }

3083: /*@C
3084:    MatLUFactorSymbolic - Performs symbolic LU factorization of matrix.
3085:    Call this routine before calling MatLUFactorNumeric().

3087:    Collective on Mat

3089:    Input Parameters:
3090: +  fact - the factor matrix obtained with MatGetFactor()
3091: .  mat - the matrix
3092: .  row, col - row and column permutations
3093: -  info - options for factorization, includes
3094: $          fill - expected fill as ratio of original fill.
3095: $          dtcol - pivot tolerance (0 no pivot, 1 full column pivoting)
3096: $                   Run with the option -info to determine an optimal value to use


3099:    Notes:
3100:     See Users-Manual: ch_mat for additional information about choosing the fill factor for better efficiency.

3102:    Most users should employ the simplified KSP interface for linear solvers
3103:    instead of working directly with matrix algebra routines such as this.
3104:    See, e.g., KSPCreate().

3106:    Level: developer

3108: .seealso: MatLUFactor(), MatLUFactorNumeric(), MatCholeskyFactor(), MatFactorInfo, MatFactorInfoInitialize()

3110:     Developer Note: fortran interface is not autogenerated as the f90
3111:     interface defintion cannot be generated correctly [due to MatFactorInfo]

3113: @*/
3114: PetscErrorCode MatLUFactorSymbolic(Mat fact,Mat mat,IS row,IS col,const MatFactorInfo *info)
3115: {
3117:   MatFactorInfo  tinfo;

3126:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
3127:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
3128:   if (!(fact)->ops->lufactorsymbolic) {
3129:     MatSolverType stype;
3130:     MatFactorGetSolverType(fact,&stype);
3131:     SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Matrix type %s symbolic LU using solver package %s",((PetscObject)mat)->type_name,stype);
3132:   }
3133:   MatCheckPreallocated(mat,2);
3134:   if (!info) {
3135:     MatFactorInfoInitialize(&tinfo);
3136:     info = &tinfo;
3137:   }

3139:   PetscLogEventBegin(MAT_LUFactorSymbolic,mat,row,col,0);
3140:   (fact->ops->lufactorsymbolic)(fact,mat,row,col,info);
3141:   PetscLogEventEnd(MAT_LUFactorSymbolic,mat,row,col,0);
3142:   PetscObjectStateIncrease((PetscObject)fact);
3143:   return(0);
3144: }

3146: /*@C
3147:    MatLUFactorNumeric - Performs numeric LU factorization of a matrix.
3148:    Call this routine after first calling MatLUFactorSymbolic().

3150:    Collective on Mat

3152:    Input Parameters:
3153: +  fact - the factor matrix obtained with MatGetFactor()
3154: .  mat - the matrix
3155: -  info - options for factorization

3157:    Notes:
3158:    See MatLUFactor() for in-place factorization.  See
3159:    MatCholeskyFactorNumeric() for the symmetric, positive definite case.

3161:    Most users should employ the simplified KSP interface for linear solvers
3162:    instead of working directly with matrix algebra routines such as this.
3163:    See, e.g., KSPCreate().

3165:    Level: developer

3167: .seealso: MatLUFactorSymbolic(), MatLUFactor(), MatCholeskyFactor()

3169:     Developer Note: fortran interface is not autogenerated as the f90
3170:     interface defintion cannot be generated correctly [due to MatFactorInfo]

3172: @*/
3173: PetscErrorCode MatLUFactorNumeric(Mat fact,Mat mat,const MatFactorInfo *info)
3174: {
3175:   MatFactorInfo  tinfo;

3183:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
3184:   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);

3186:   if (!(fact)->ops->lufactornumeric) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s numeric LU",((PetscObject)mat)->type_name);
3187:   MatCheckPreallocated(mat,2);
3188:   if (!info) {
3189:     MatFactorInfoInitialize(&tinfo);
3190:     info = &tinfo;
3191:   }

3193:   PetscLogEventBegin(MAT_LUFactorNumeric,mat,fact,0,0);
3194:   (fact->ops->lufactornumeric)(fact,mat,info);
3195:   PetscLogEventEnd(MAT_LUFactorNumeric,mat,fact,0,0);
3196:   MatViewFromOptions(fact,NULL,"-mat_factor_view");
3197:   PetscObjectStateIncrease((PetscObject)fact);
3198:   return(0);
3199: }

3201: /*@C
3202:    MatCholeskyFactor - Performs in-place Cholesky factorization of a
3203:    symmetric matrix.

3205:    Collective on Mat

3207:    Input Parameters:
3208: +  mat - the matrix
3209: .  perm - row and column permutations
3210: -  f - expected fill as ratio of original fill

3212:    Notes:
3213:    See MatLUFactor() for the nonsymmetric case.  See also
3214:    MatCholeskyFactorSymbolic(), and MatCholeskyFactorNumeric().

3216:    Most users should employ the simplified KSP interface for linear solvers
3217:    instead of working directly with matrix algebra routines such as this.
3218:    See, e.g., KSPCreate().

3220:    Level: developer

3222: .seealso: MatLUFactor(), MatCholeskyFactorSymbolic(), MatCholeskyFactorNumeric()
3223:           MatGetOrdering()

3225:     Developer Note: fortran interface is not autogenerated as the f90
3226:     interface defintion cannot be generated correctly [due to MatFactorInfo]

3228: @*/
3229: PetscErrorCode MatCholeskyFactor(Mat mat,IS perm,const MatFactorInfo *info)
3230: {
3232:   MatFactorInfo  tinfo;

3239:   if (mat->rmap->N != mat->cmap->N) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONG,"Matrix must be square");
3240:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
3241:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
3242:   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);
3243:   MatCheckPreallocated(mat,1);
3244:   if (!info) {
3245:     MatFactorInfoInitialize(&tinfo);
3246:     info = &tinfo;
3247:   }

3249:   PetscLogEventBegin(MAT_CholeskyFactor,mat,perm,0,0);
3250:   (*mat->ops->choleskyfactor)(mat,perm,info);
3251:   PetscLogEventEnd(MAT_CholeskyFactor,mat,perm,0,0);
3252:   PetscObjectStateIncrease((PetscObject)mat);
3253:   return(0);
3254: }

3256: /*@C
3257:    MatCholeskyFactorSymbolic - Performs symbolic Cholesky factorization
3258:    of a symmetric matrix.

3260:    Collective on Mat

3262:    Input Parameters:
3263: +  fact - the factor matrix obtained with MatGetFactor()
3264: .  mat - the matrix
3265: .  perm - row and column permutations
3266: -  info - options for factorization, includes
3267: $          fill - expected fill as ratio of original fill.
3268: $          dtcol - pivot tolerance (0 no pivot, 1 full column pivoting)
3269: $                   Run with the option -info to determine an optimal value to use

3271:    Notes:
3272:    See MatLUFactorSymbolic() for the nonsymmetric case.  See also
3273:    MatCholeskyFactor() and MatCholeskyFactorNumeric().

3275:    Most users should employ the simplified KSP interface for linear solvers
3276:    instead of working directly with matrix algebra routines such as this.
3277:    See, e.g., KSPCreate().

3279:    Level: developer

3281: .seealso: MatLUFactorSymbolic(), MatCholeskyFactor(), MatCholeskyFactorNumeric()
3282:           MatGetOrdering()

3284:     Developer Note: fortran interface is not autogenerated as the f90
3285:     interface defintion cannot be generated correctly [due to MatFactorInfo]

3287: @*/
3288: PetscErrorCode MatCholeskyFactorSymbolic(Mat fact,Mat mat,IS perm,const MatFactorInfo *info)
3289: {
3291:   MatFactorInfo  tinfo;

3299:   if (mat->rmap->N != mat->cmap->N) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONG,"Matrix must be square");
3300:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
3301:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
3302:   if (!(fact)->ops->choleskyfactorsymbolic) {
3303:     MatSolverType stype;
3304:     MatFactorGetSolverType(fact,&stype);
3305:     SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s symbolic factor Cholesky using solver package %s",((PetscObject)mat)->type_name,stype);
3306:   }
3307:   MatCheckPreallocated(mat,2);
3308:   if (!info) {
3309:     MatFactorInfoInitialize(&tinfo);
3310:     info = &tinfo;
3311:   }

3313:   PetscLogEventBegin(MAT_CholeskyFactorSymbolic,mat,perm,0,0);
3314:   (fact->ops->choleskyfactorsymbolic)(fact,mat,perm,info);
3315:   PetscLogEventEnd(MAT_CholeskyFactorSymbolic,mat,perm,0,0);
3316:   PetscObjectStateIncrease((PetscObject)fact);
3317:   return(0);
3318: }

3320: /*@C
3321:    MatCholeskyFactorNumeric - Performs numeric Cholesky factorization
3322:    of a symmetric matrix. Call this routine after first calling
3323:    MatCholeskyFactorSymbolic().

3325:    Collective on Mat

3327:    Input Parameters:
3328: +  fact - the factor matrix obtained with MatGetFactor()
3329: .  mat - the initial matrix
3330: .  info - options for factorization
3331: -  fact - the symbolic factor of mat


3334:    Notes:
3335:    Most users should employ the simplified KSP interface for linear solvers
3336:    instead of working directly with matrix algebra routines such as this.
3337:    See, e.g., KSPCreate().

3339:    Level: developer

3341: .seealso: MatCholeskyFactorSymbolic(), MatCholeskyFactor(), MatLUFactorNumeric()

3343:     Developer Note: fortran interface is not autogenerated as the f90
3344:     interface defintion cannot be generated correctly [due to MatFactorInfo]

3346: @*/
3347: PetscErrorCode MatCholeskyFactorNumeric(Mat fact,Mat mat,const MatFactorInfo *info)
3348: {
3349:   MatFactorInfo  tinfo;

3357:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
3358:   if (!(fact)->ops->choleskyfactornumeric) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s numeric factor Cholesky",((PetscObject)mat)->type_name);
3359:   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);
3360:   MatCheckPreallocated(mat,2);
3361:   if (!info) {
3362:     MatFactorInfoInitialize(&tinfo);
3363:     info = &tinfo;
3364:   }

3366:   PetscLogEventBegin(MAT_CholeskyFactorNumeric,mat,fact,0,0);
3367:   (fact->ops->choleskyfactornumeric)(fact,mat,info);
3368:   PetscLogEventEnd(MAT_CholeskyFactorNumeric,mat,fact,0,0);
3369:   MatViewFromOptions(fact,NULL,"-mat_factor_view");
3370:   PetscObjectStateIncrease((PetscObject)fact);
3371:   return(0);
3372: }

3374: /*@
3375:    MatQRFactor - Performs in-place QR factorization of matrix.

3377:    Collective on Mat

3379:    Input Parameters:
3380: +  mat - the matrix
3381: .  col - column permutation
3382: -  info - options for factorization, includes
3383: $          fill - expected fill as ratio of original fill.
3384: $          dtcol - pivot tolerance (0 no pivot, 1 full column pivoting)
3385: $                   Run with the option -info to determine an optimal value to use

3387:    Notes:
3388:    Most users should employ the simplified KSP interface for linear solvers
3389:    instead of working directly with matrix algebra routines such as this.
3390:    See, e.g., KSPCreate().

3392:    This changes the state of the matrix to a factored matrix; it cannot be used
3393:    for example with MatSetValues() unless one first calls MatSetUnfactored().

3395:    Level: developer

3397: .seealso: MatQRFactorSymbolic(), MatQRFactorNumeric(), MatLUFactor(),
3398:           MatSetUnfactored(), MatFactorInfo, MatGetFactor()

3400:     Developer Note: fortran interface is not autogenerated as the f90
3401:     interface defintion cannot be generated correctly [due to MatFactorInfo]

3403: @*/
3404: PetscErrorCode MatQRFactor(Mat mat, IS col, const MatFactorInfo *info)
3405: {

3413:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
3414:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
3415:   MatCheckPreallocated(mat,1);
3416:   PetscLogEventBegin(MAT_QRFactor,mat,col,0,0);
3417:   PetscUseMethod(mat,"MatQRFactor_C", (Mat,IS,const MatFactorInfo*), (mat, col, info));
3418:   PetscLogEventEnd(MAT_QRFactor,mat,col,0,0);
3419:   PetscObjectStateIncrease((PetscObject)mat);
3420:   return(0);
3421: }

3423: /*@
3424:    MatQRFactorSymbolic - Performs symbolic QR factorization of matrix.
3425:    Call this routine before calling MatQRFactorNumeric().

3427:    Collective on Mat

3429:    Input Parameters:
3430: +  fact - the factor matrix obtained with MatGetFactor()
3431: .  mat - the matrix
3432: .  col - column permutation
3433: -  info - options for factorization, includes
3434: $          fill - expected fill as ratio of original fill.
3435: $          dtcol - pivot tolerance (0 no pivot, 1 full column pivoting)
3436: $                   Run with the option -info to determine an optimal value to use

3438:    Most users should employ the simplified KSP interface for linear solvers
3439:    instead of working directly with matrix algebra routines such as this.
3440:    See, e.g., KSPCreate().

3442:    Level: developer

3444: .seealso: MatQRFactor(), MatQRFactorNumeric(), MatLUFactor(), MatFactorInfo, MatFactorInfoInitialize()

3446:     Developer Note: fortran interface is not autogenerated as the f90
3447:     interface defintion cannot be generated correctly [due to MatFactorInfo]

3449: @*/
3450: PetscErrorCode MatQRFactorSymbolic(Mat fact,Mat mat,IS col,const MatFactorInfo *info)
3451: {
3453:   MatFactorInfo  tinfo;

3461:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
3462:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
3463:   MatCheckPreallocated(mat,2);
3464:   if (!info) {
3465:     MatFactorInfoInitialize(&tinfo);
3466:     info = &tinfo;
3467:   }

3469:   PetscLogEventBegin(MAT_QRFactorSymbolic,fact,mat,col,0);
3470:   PetscUseMethod(fact,"MatQRFactorSymbolic_C", (Mat,Mat,IS,const MatFactorInfo*), (fact, mat, col, info));
3471:   PetscLogEventEnd(MAT_QRFactorSymbolic,fact,mat,col,0);
3472:   PetscObjectStateIncrease((PetscObject)fact);
3473:   return(0);
3474: }

3476: /*@
3477:    MatQRFactorNumeric - Performs numeric QR factorization of a matrix.
3478:    Call this routine after first calling MatQRFactorSymbolic().

3480:    Collective on Mat

3482:    Input Parameters:
3483: +  fact - the factor matrix obtained with MatGetFactor()
3484: .  mat - the matrix
3485: -  info - options for factorization

3487:    Notes:
3488:    See MatQRFactor() for in-place factorization.

3490:    Most users should employ the simplified KSP interface for linear solvers
3491:    instead of working directly with matrix algebra routines such as this.
3492:    See, e.g., KSPCreate().

3494:    Level: developer

3496: .seealso: MatQRFactorSymbolic(), MatLUFactor()

3498:     Developer Note: fortran interface is not autogenerated as the f90
3499:     interface defintion cannot be generated correctly [due to MatFactorInfo]

3501: @*/
3502: PetscErrorCode MatQRFactorNumeric(Mat fact,Mat mat,const MatFactorInfo *info)
3503: {
3504:   MatFactorInfo  tinfo;

3512:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
3513:   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);

3515:   MatCheckPreallocated(mat,2);
3516:   if (!info) {
3517:     MatFactorInfoInitialize(&tinfo);
3518:     info = &tinfo;
3519:   }

3521:   PetscLogEventBegin(MAT_QRFactorNumeric,mat,fact,0,0);
3522:   PetscUseMethod(fact,"MatQRFactorNumeric_C", (Mat,Mat,const MatFactorInfo*), (fact, mat, info));
3523:   PetscLogEventEnd(MAT_QRFactorNumeric,mat,fact,0,0);
3524:   MatViewFromOptions(fact,NULL,"-mat_factor_view");
3525:   PetscObjectStateIncrease((PetscObject)fact);
3526:   return(0);
3527: }

3529: /* ----------------------------------------------------------------*/
3530: /*@
3531:    MatSolve - Solves A x = b, given a factored matrix.

3533:    Neighbor-wise Collective on Mat

3535:    Input Parameters:
3536: +  mat - the factored matrix
3537: -  b - the right-hand-side vector

3539:    Output Parameter:
3540: .  x - the result vector

3542:    Notes:
3543:    The vectors b and x cannot be the same.  I.e., one cannot
3544:    call MatSolve(A,x,x).

3546:    Notes:
3547:    Most users should employ the simplified KSP interface for linear solvers
3548:    instead of working directly with matrix algebra routines such as this.
3549:    See, e.g., KSPCreate().

3551:    Level: developer

3553: .seealso: MatSolveAdd(), MatSolveTranspose(), MatSolveTransposeAdd()
3554: @*/
3555: PetscErrorCode MatSolve(Mat mat,Vec b,Vec x)
3556: {

3566:   if (x == b) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_IDN,"x and b must be different vectors");
3567:   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);
3568:   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);
3569:   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);
3570:   if (!mat->rmap->N && !mat->cmap->N) return(0);
3571:   MatCheckPreallocated(mat,1);

3573:   PetscLogEventBegin(MAT_Solve,mat,b,x,0);
3574:   if (mat->factorerrortype) {
3575:     PetscInfo1(mat,"MatFactorError %D\n",mat->factorerrortype);
3576:     VecSetInf(x);
3577:   } else {
3578:     if (!mat->ops->solve) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
3579:     (*mat->ops->solve)(mat,b,x);
3580:   }
3581:   PetscLogEventEnd(MAT_Solve,mat,b,x,0);
3582:   PetscObjectStateIncrease((PetscObject)x);
3583:   return(0);
3584: }

3586: static PetscErrorCode MatMatSolve_Basic(Mat A,Mat B,Mat X,PetscBool trans)
3587: {
3589:   Vec            b,x;
3590:   PetscInt       m,N,i;
3591:   PetscScalar    *bb,*xx;
3592:   PetscErrorCode (*f)(Mat,Vec,Vec);

3595:   if (A->factorerrortype) {
3596:     PetscInfo1(A,"MatFactorError %D\n",A->factorerrortype);
3597:     MatSetInf(X);
3598:     return(0);
3599:   }
3600:   f = trans ? A->ops->solvetranspose : A->ops->solve;
3601:   if (!f) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Mat type %s",((PetscObject)A)->type_name);

3603:   MatDenseGetArrayRead(B,(const PetscScalar**)&bb);
3604:   MatDenseGetArray(X,&xx);
3605:   MatGetLocalSize(B,&m,NULL);  /* number local rows */
3606:   MatGetSize(B,NULL,&N);       /* total columns in dense matrix */
3607:   MatCreateVecs(A,&x,&b);
3608:   for (i=0; i<N; i++) {
3609:     VecPlaceArray(b,bb + i*m);
3610:     VecPlaceArray(x,xx + i*m);
3611:     (*f)(A,b,x);
3612:     VecResetArray(x);
3613:     VecResetArray(b);
3614:   }
3615:   VecDestroy(&b);
3616:   VecDestroy(&x);
3617:   MatDenseRestoreArrayRead(B,(const PetscScalar**)&bb);
3618:   MatDenseRestoreArray(X,&xx);
3619:   return(0);
3620: }

3622: /*@
3623:    MatMatSolve - Solves A X = B, given a factored matrix.

3625:    Neighbor-wise Collective on Mat

3627:    Input Parameters:
3628: +  A - the factored matrix
3629: -  B - the right-hand-side matrix MATDENSE (or sparse -- when using MUMPS)

3631:    Output Parameter:
3632: .  X - the result matrix (dense matrix)

3634:    Notes:
3635:    If B is a MATDENSE matrix then one can call MatMatSolve(A,B,B) except with MKL_CPARDISO;
3636:    otherwise, B and X cannot be the same.

3638:    Notes:
3639:    Most users should usually employ the simplified KSP interface for linear solvers
3640:    instead of working directly with matrix algebra routines such as this.
3641:    See, e.g., KSPCreate(). However KSP can only solve for one vector (column of X)
3642:    at a time.

3644:    Level: developer

3646: .seealso: MatMatSolveTranspose(), MatLUFactor(), MatCholeskyFactor()
3647: @*/
3648: PetscErrorCode MatMatSolve(Mat A,Mat B,Mat X)
3649: {

3659:   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);
3660:   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);
3661:   if (X->cmap->N != B->cmap->N) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Solution matrix must have same number of columns as rhs matrix");
3662:   if (!A->rmap->N && !A->cmap->N) return(0);
3663:   if (!A->factortype) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Unfactored matrix");
3664:   MatCheckPreallocated(A,1);

3666:   PetscLogEventBegin(MAT_MatSolve,A,B,X,0);
3667:   if (!A->ops->matsolve) {
3668:     PetscInfo1(A,"Mat type %s using basic MatMatSolve\n",((PetscObject)A)->type_name);
3669:     MatMatSolve_Basic(A,B,X,PETSC_FALSE);
3670:   } else {
3671:     (*A->ops->matsolve)(A,B,X);
3672:   }
3673:   PetscLogEventEnd(MAT_MatSolve,A,B,X,0);
3674:   PetscObjectStateIncrease((PetscObject)X);
3675:   return(0);
3676: }

3678: /*@
3679:    MatMatSolveTranspose - Solves A^T X = B, given a factored matrix.

3681:    Neighbor-wise Collective on Mat

3683:    Input Parameters:
3684: +  A - the factored matrix
3685: -  B - the right-hand-side matrix  (dense matrix)

3687:    Output Parameter:
3688: .  X - the result matrix (dense matrix)

3690:    Notes:
3691:    The matrices B and X cannot be the same.  I.e., one cannot
3692:    call MatMatSolveTranspose(A,X,X).

3694:    Notes:
3695:    Most users should usually employ the simplified KSP interface for linear solvers
3696:    instead of working directly with matrix algebra routines such as this.
3697:    See, e.g., KSPCreate(). However KSP can only solve for one vector (column of X)
3698:    at a time.

3700:    When using SuperLU_Dist or MUMPS as a parallel solver, PETSc will use their functionality to solve multiple right hand sides simultaneously.

3702:    Level: developer

3704: .seealso: MatMatSolve(), MatLUFactor(), MatCholeskyFactor()
3705: @*/
3706: PetscErrorCode MatMatSolveTranspose(Mat A,Mat B,Mat X)
3707: {

3717:   if (X == B) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_IDN,"X and B must be different matrices");
3718:   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);
3719:   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);
3720:   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);
3721:   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");
3722:   if (!A->rmap->N && !A->cmap->N) return(0);
3723:   if (!A->factortype) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Unfactored matrix");
3724:   MatCheckPreallocated(A,1);

3726:   PetscLogEventBegin(MAT_MatSolve,A,B,X,0);
3727:   if (!A->ops->matsolvetranspose) {
3728:     PetscInfo1(A,"Mat type %s using basic MatMatSolveTranspose\n",((PetscObject)A)->type_name);
3729:     MatMatSolve_Basic(A,B,X,PETSC_TRUE);
3730:   } else {
3731:     (*A->ops->matsolvetranspose)(A,B,X);
3732:   }
3733:   PetscLogEventEnd(MAT_MatSolve,A,B,X,0);
3734:   PetscObjectStateIncrease((PetscObject)X);
3735:   return(0);
3736: }

3738: /*@
3739:    MatMatTransposeSolve - Solves A X = B^T, given a factored matrix.

3741:    Neighbor-wise Collective on Mat

3743:    Input Parameters:
3744: +  A - the factored matrix
3745: -  Bt - the transpose of right-hand-side matrix

3747:    Output Parameter:
3748: .  X - the result matrix (dense matrix)

3750:    Notes:
3751:    Most users should usually employ the simplified KSP interface for linear solvers
3752:    instead of working directly with matrix algebra routines such as this.
3753:    See, e.g., KSPCreate(). However KSP can only solve for one vector (column of X)
3754:    at a time.

3756:    For MUMPS, it only supports centralized sparse compressed column format on the host processor for right hand side matrix. User must create B^T in sparse compressed row format on the host processor and call MatMatTransposeSolve() to implement MUMPS' MatMatSolve().

3758:    Level: developer

3760: .seealso: MatMatSolve(), MatMatSolveTranspose(), MatLUFactor(), MatCholeskyFactor()
3761: @*/
3762: PetscErrorCode MatMatTransposeSolve(Mat A,Mat Bt,Mat X)
3763: {


3774:   if (X == Bt) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_IDN,"X and B must be different matrices");
3775:   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);
3776:   if (A->rmap->N != Bt->cmap->N) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Mat A,Mat Bt: global dim %D %D",A->rmap->N,Bt->cmap->N);
3777:   if (X->cmap->N < Bt->rmap->N) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_SIZ,"Solution matrix must have same number of columns as row number of the rhs matrix");
3778:   if (!A->rmap->N && !A->cmap->N) return(0);
3779:   if (!A->factortype) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Unfactored matrix");
3780:   MatCheckPreallocated(A,1);

3782:   if (!A->ops->mattransposesolve) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Mat type %s",((PetscObject)A)->type_name);
3783:   PetscLogEventBegin(MAT_MatTrSolve,A,Bt,X,0);
3784:   (*A->ops->mattransposesolve)(A,Bt,X);
3785:   PetscLogEventEnd(MAT_MatTrSolve,A,Bt,X,0);
3786:   PetscObjectStateIncrease((PetscObject)X);
3787:   return(0);
3788: }

3790: /*@
3791:    MatForwardSolve - Solves L x = b, given a factored matrix, A = LU, or
3792:                             U^T*D^(1/2) x = b, given a factored symmetric matrix, A = U^T*D*U,

3794:    Neighbor-wise Collective on Mat

3796:    Input Parameters:
3797: +  mat - the factored matrix
3798: -  b - the right-hand-side vector

3800:    Output Parameter:
3801: .  x - the result vector

3803:    Notes:
3804:    MatSolve() should be used for most applications, as it performs
3805:    a forward solve followed by a backward solve.

3807:    The vectors b and x cannot be the same,  i.e., one cannot
3808:    call MatForwardSolve(A,x,x).

3810:    For matrix in seqsbaij format with block size larger than 1,
3811:    the diagonal blocks are not implemented as D = D^(1/2) * D^(1/2) yet.
3812:    MatForwardSolve() solves U^T*D y = b, and
3813:    MatBackwardSolve() solves U x = y.
3814:    Thus they do not provide a symmetric preconditioner.

3816:    Most users should employ the simplified KSP interface for linear solvers
3817:    instead of working directly with matrix algebra routines such as this.
3818:    See, e.g., KSPCreate().

3820:    Level: developer

3822: .seealso: MatSolve(), MatBackwardSolve()
3823: @*/
3824: PetscErrorCode MatForwardSolve(Mat mat,Vec b,Vec x)
3825: {

3835:   if (x == b) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_IDN,"x and b must be different vectors");
3836:   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);
3837:   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);
3838:   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);
3839:   if (!mat->rmap->N && !mat->cmap->N) return(0);
3840:   MatCheckPreallocated(mat,1);

3842:   if (!mat->ops->forwardsolve) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
3843:   PetscLogEventBegin(MAT_ForwardSolve,mat,b,x,0);
3844:   (*mat->ops->forwardsolve)(mat,b,x);
3845:   PetscLogEventEnd(MAT_ForwardSolve,mat,b,x,0);
3846:   PetscObjectStateIncrease((PetscObject)x);
3847:   return(0);
3848: }

3850: /*@
3851:    MatBackwardSolve - Solves U x = b, given a factored matrix, A = LU.
3852:                              D^(1/2) U x = b, given a factored symmetric matrix, A = U^T*D*U,

3854:    Neighbor-wise Collective on Mat

3856:    Input Parameters:
3857: +  mat - the factored matrix
3858: -  b - the right-hand-side vector

3860:    Output Parameter:
3861: .  x - the result vector

3863:    Notes:
3864:    MatSolve() should be used for most applications, as it performs
3865:    a forward solve followed by a backward solve.

3867:    The vectors b and x cannot be the same.  I.e., one cannot
3868:    call MatBackwardSolve(A,x,x).

3870:    For matrix in seqsbaij format with block size larger than 1,
3871:    the diagonal blocks are not implemented as D = D^(1/2) * D^(1/2) yet.
3872:    MatForwardSolve() solves U^T*D y = b, and
3873:    MatBackwardSolve() solves U x = y.
3874:    Thus they do not provide a symmetric preconditioner.

3876:    Most users should employ the simplified KSP interface for linear solvers
3877:    instead of working directly with matrix algebra routines such as this.
3878:    See, e.g., KSPCreate().

3880:    Level: developer

3882: .seealso: MatSolve(), MatForwardSolve()
3883: @*/
3884: PetscErrorCode MatBackwardSolve(Mat mat,Vec b,Vec x)
3885: {

3895:   if (x == b) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_IDN,"x and b must be different vectors");
3896:   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);
3897:   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);
3898:   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);
3899:   if (!mat->rmap->N && !mat->cmap->N) return(0);
3900:   MatCheckPreallocated(mat,1);

3902:   if (!mat->ops->backwardsolve) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
3903:   PetscLogEventBegin(MAT_BackwardSolve,mat,b,x,0);
3904:   (*mat->ops->backwardsolve)(mat,b,x);
3905:   PetscLogEventEnd(MAT_BackwardSolve,mat,b,x,0);
3906:   PetscObjectStateIncrease((PetscObject)x);
3907:   return(0);
3908: }

3910: /*@
3911:    MatSolveAdd - Computes x = y + inv(A)*b, given a factored matrix.

3913:    Neighbor-wise Collective on Mat

3915:    Input Parameters:
3916: +  mat - the factored matrix
3917: .  b - the right-hand-side vector
3918: -  y - the vector to be added to

3920:    Output Parameter:
3921: .  x - the result vector

3923:    Notes:
3924:    The vectors b and x cannot be the same.  I.e., one cannot
3925:    call MatSolveAdd(A,x,y,x).

3927:    Most users should employ the simplified KSP interface for linear solvers
3928:    instead of working directly with matrix algebra routines such as this.
3929:    See, e.g., KSPCreate().

3931:    Level: developer

3933: .seealso: MatSolve(), MatSolveTranspose(), MatSolveTransposeAdd()
3934: @*/
3935: PetscErrorCode MatSolveAdd(Mat mat,Vec b,Vec y,Vec x)
3936: {
3937:   PetscScalar    one = 1.0;
3938:   Vec            tmp;

3950:   if (x == b) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_IDN,"x and b must be different vectors");
3951:   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);
3952:   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);
3953:   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);
3954:   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);
3955:   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);
3956:   if (!mat->rmap->N && !mat->cmap->N) return(0);
3957:    MatCheckPreallocated(mat,1);

3959:   PetscLogEventBegin(MAT_SolveAdd,mat,b,x,y);
3960:   if (mat->factorerrortype) {
3961:     PetscInfo1(mat,"MatFactorError %D\n",mat->factorerrortype);
3962:     VecSetInf(x);
3963:   } else if (mat->ops->solveadd) {
3964:     (*mat->ops->solveadd)(mat,b,y,x);
3965:   } else {
3966:     /* do the solve then the add manually */
3967:     if (x != y) {
3968:       MatSolve(mat,b,x);
3969:       VecAXPY(x,one,y);
3970:     } else {
3971:       VecDuplicate(x,&tmp);
3972:       PetscLogObjectParent((PetscObject)mat,(PetscObject)tmp);
3973:       VecCopy(x,tmp);
3974:       MatSolve(mat,b,x);
3975:       VecAXPY(x,one,tmp);
3976:       VecDestroy(&tmp);
3977:     }
3978:   }
3979:   PetscLogEventEnd(MAT_SolveAdd,mat,b,x,y);
3980:   PetscObjectStateIncrease((PetscObject)x);
3981:   return(0);
3982: }

3984: /*@
3985:    MatSolveTranspose - Solves A' x = b, given a factored matrix.

3987:    Neighbor-wise Collective on Mat

3989:    Input Parameters:
3990: +  mat - the factored matrix
3991: -  b - the right-hand-side vector

3993:    Output Parameter:
3994: .  x - the result vector

3996:    Notes:
3997:    The vectors b and x cannot be the same.  I.e., one cannot
3998:    call MatSolveTranspose(A,x,x).

4000:    Most users should employ the simplified KSP interface for linear solvers
4001:    instead of working directly with matrix algebra routines such as this.
4002:    See, e.g., KSPCreate().

4004:    Level: developer

4006: .seealso: MatSolve(), MatSolveAdd(), MatSolveTransposeAdd()
4007: @*/
4008: PetscErrorCode MatSolveTranspose(Mat mat,Vec b,Vec x)
4009: {

4019:   if (x == b) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_IDN,"x and b must be different vectors");
4020:   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);
4021:   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);
4022:   if (!mat->rmap->N && !mat->cmap->N) return(0);
4023:   MatCheckPreallocated(mat,1);
4024:   PetscLogEventBegin(MAT_SolveTranspose,mat,b,x,0);
4025:   if (mat->factorerrortype) {
4026:     PetscInfo1(mat,"MatFactorError %D\n",mat->factorerrortype);
4027:     VecSetInf(x);
4028:   } else {
4029:     if (!mat->ops->solvetranspose) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Matrix type %s",((PetscObject)mat)->type_name);
4030:     (*mat->ops->solvetranspose)(mat,b,x);
4031:   }
4032:   PetscLogEventEnd(MAT_SolveTranspose,mat,b,x,0);
4033:   PetscObjectStateIncrease((PetscObject)x);
4034:   return(0);
4035: }

4037: /*@
4038:    MatSolveTransposeAdd - Computes x = y + inv(Transpose(A)) b, given a
4039:                       factored matrix.

4041:    Neighbor-wise Collective on Mat

4043:    Input Parameters:
4044: +  mat - the factored matrix
4045: .  b - the right-hand-side vector
4046: -  y - the vector to be added to

4048:    Output Parameter:
4049: .  x - the result vector

4051:    Notes:
4052:    The vectors b and x cannot be the same.  I.e., one cannot
4053:    call MatSolveTransposeAdd(A,x,y,x).

4055:    Most users should employ the simplified KSP interface for linear solvers
4056:    instead of working directly with matrix algebra routines such as this.
4057:    See, e.g., KSPCreate().

4059:    Level: developer

4061: .seealso: MatSolve(), MatSolveAdd(), MatSolveTranspose()
4062: @*/
4063: PetscErrorCode MatSolveTransposeAdd(Mat mat,Vec b,Vec y,Vec x)
4064: {
4065:   PetscScalar    one = 1.0;
4067:   Vec            tmp;

4078:   if (x == b) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_IDN,"x and b must be different vectors");
4079:   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);
4080:   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);
4081:   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);
4082:   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);
4083:   if (!mat->rmap->N && !mat->cmap->N) return(0);
4084:    MatCheckPreallocated(mat,1);

4086:   PetscLogEventBegin(MAT_SolveTransposeAdd,mat,b,x,y);
4087:   if (mat->factorerrortype) {
4088:     PetscInfo1(mat,"MatFactorError %D\n",mat->factorerrortype);
4089:     VecSetInf(x);
4090:   } else if (mat->ops->solvetransposeadd){
4091:     (*mat->ops->solvetransposeadd)(mat,b,y,x);
4092:   } else {
4093:     /* do the solve then the add manually */
4094:     if (x != y) {
4095:       MatSolveTranspose(mat,b,x);
4096:       VecAXPY(x,one,y);
4097:     } else {
4098:       VecDuplicate(x,&tmp);
4099:       PetscLogObjectParent((PetscObject)mat,(PetscObject)tmp);
4100:       VecCopy(x,tmp);
4101:       MatSolveTranspose(mat,b,x);
4102:       VecAXPY(x,one,tmp);
4103:       VecDestroy(&tmp);
4104:     }
4105:   }
4106:   PetscLogEventEnd(MAT_SolveTransposeAdd,mat,b,x,y);
4107:   PetscObjectStateIncrease((PetscObject)x);
4108:   return(0);
4109: }
4110: /* ----------------------------------------------------------------*/

4112: /*@
4113:    MatSOR - Computes relaxation (SOR, Gauss-Seidel) sweeps.

4115:    Neighbor-wise Collective on Mat

4117:    Input Parameters:
4118: +  mat - the matrix
4119: .  b - the right hand side
4120: .  omega - the relaxation factor
4121: .  flag - flag indicating the type of SOR (see below)
4122: .  shift -  diagonal shift
4123: .  its - the number of iterations
4124: -  lits - the number of local iterations

4126:    Output Parameters:
4127: .  x - the solution (can contain an initial guess, use option SOR_ZERO_INITIAL_GUESS to indicate no guess)

4129:    SOR Flags:
4130: +     SOR_FORWARD_SWEEP - forward SOR
4131: .     SOR_BACKWARD_SWEEP - backward SOR
4132: .     SOR_SYMMETRIC_SWEEP - SSOR (symmetric SOR)
4133: .     SOR_LOCAL_FORWARD_SWEEP - local forward SOR
4134: .     SOR_LOCAL_BACKWARD_SWEEP - local forward SOR
4135: .     SOR_LOCAL_SYMMETRIC_SWEEP - local SSOR
4136: .     SOR_APPLY_UPPER, SOR_APPLY_LOWER - applies
4137:          upper/lower triangular part of matrix to
4138:          vector (with omega)
4139: -     SOR_ZERO_INITIAL_GUESS - zero initial guess

4141:    Notes:
4142:    SOR_LOCAL_FORWARD_SWEEP, SOR_LOCAL_BACKWARD_SWEEP, and
4143:    SOR_LOCAL_SYMMETRIC_SWEEP perform separate independent smoothings
4144:    on each processor.

4146:    Application programmers will not generally use MatSOR() directly,
4147:    but instead will employ the KSP/PC interface.

4149:    Notes:
4150:     for BAIJ, SBAIJ, and AIJ matrices with Inodes this does a block SOR smoothing, otherwise it does a pointwise smoothing

4152:    Notes for Advanced Users:
4153:    The flags are implemented as bitwise inclusive or operations.
4154:    For example, use (SOR_ZERO_INITIAL_GUESS | SOR_SYMMETRIC_SWEEP)
4155:    to specify a zero initial guess for SSOR.

4157:    Most users should employ the simplified KSP interface for linear solvers
4158:    instead of working directly with matrix algebra routines such as this.
4159:    See, e.g., KSPCreate().

4161:    Vectors x and b CANNOT be the same

4163:    Developer Note: We should add block SOR support for AIJ matrices with block size set to great than one and no inodes

4165:    Level: developer

4167: @*/
4168: PetscErrorCode MatSOR(Mat mat,Vec b,PetscReal omega,MatSORType flag,PetscReal shift,PetscInt its,PetscInt lits,Vec x)
4169: {

4179:   if (!mat->ops->sor) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
4180:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
4181:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
4182:   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);
4183:   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);
4184:   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);
4185:   if (its <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D positive",its);
4186:   if (lits <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires local its %D positive",lits);
4187:   if (b == x) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_IDN,"b and x vector cannot be the same");

4189:   MatCheckPreallocated(mat,1);
4190:   PetscLogEventBegin(MAT_SOR,mat,b,x,0);
4191:   ierr =(*mat->ops->sor)(mat,b,omega,flag,shift,its,lits,x);
4192:   PetscLogEventEnd(MAT_SOR,mat,b,x,0);
4193:   PetscObjectStateIncrease((PetscObject)x);
4194:   return(0);
4195: }

4197: /*
4198:       Default matrix copy routine.
4199: */
4200: PetscErrorCode MatCopy_Basic(Mat A,Mat B,MatStructure str)
4201: {
4202:   PetscErrorCode    ierr;
4203:   PetscInt          i,rstart = 0,rend = 0,nz;
4204:   const PetscInt    *cwork;
4205:   const PetscScalar *vwork;

4208:   if (B->assembled) {
4209:     MatZeroEntries(B);
4210:   }
4211:   if (str == SAME_NONZERO_PATTERN) {
4212:     MatGetOwnershipRange(A,&rstart,&rend);
4213:     for (i=rstart; i<rend; i++) {
4214:       MatGetRow(A,i,&nz,&cwork,&vwork);
4215:       MatSetValues(B,1,&i,nz,cwork,vwork,INSERT_VALUES);
4216:       MatRestoreRow(A,i,&nz,&cwork,&vwork);
4217:     }
4218:   } else {
4219:     MatAYPX(B,0.0,A,str);
4220:   }
4221:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
4222:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
4223:   return(0);
4224: }

4226: /*@
4227:    MatCopy - Copies a matrix to another matrix.

4229:    Collective on Mat

4231:    Input Parameters:
4232: +  A - the matrix
4233: -  str - SAME_NONZERO_PATTERN or DIFFERENT_NONZERO_PATTERN

4235:    Output Parameter:
4236: .  B - where the copy is put

4238:    Notes:
4239:    If you use SAME_NONZERO_PATTERN then the two matrices must have the same nonzero pattern or the routine will crash.

4241:    MatCopy() copies the matrix entries of a matrix to another existing
4242:    matrix (after first zeroing the second matrix).  A related routine is
4243:    MatConvert(), which first creates a new matrix and then copies the data.

4245:    Level: intermediate

4247: .seealso: MatConvert(), MatDuplicate()

4249: @*/
4250: PetscErrorCode MatCopy(Mat A,Mat B,MatStructure str)
4251: {
4253:   PetscInt       i;

4261:   MatCheckPreallocated(B,2);
4262:   if (!A->assembled) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
4263:   if (A->factortype) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
4264:   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);
4265:   MatCheckPreallocated(A,1);
4266:   if (A == B) return(0);

4268:   PetscLogEventBegin(MAT_Copy,A,B,0,0);
4269:   if (A->ops->copy) {
4270:     (*A->ops->copy)(A,B,str);
4271:   } else { /* generic conversion */
4272:     MatCopy_Basic(A,B,str);
4273:   }

4275:   B->stencil.dim = A->stencil.dim;
4276:   B->stencil.noc = A->stencil.noc;
4277:   for (i=0; i<=A->stencil.dim; i++) {
4278:     B->stencil.dims[i]   = A->stencil.dims[i];
4279:     B->stencil.starts[i] = A->stencil.starts[i];
4280:   }

4282:   PetscLogEventEnd(MAT_Copy,A,B,0,0);
4283:   PetscObjectStateIncrease((PetscObject)B);
4284:   return(0);
4285: }

4287: /*@C
4288:    MatConvert - Converts a matrix to another matrix, either of the same
4289:    or different type.

4291:    Collective on Mat

4293:    Input Parameters:
4294: +  mat - the matrix
4295: .  newtype - new matrix type.  Use MATSAME to create a new matrix of the
4296:    same type as the original matrix.
4297: -  reuse - denotes if the destination matrix is to be created or reused.
4298:    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
4299:    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).

4301:    Output Parameter:
4302: .  M - pointer to place new matrix

4304:    Notes:
4305:    MatConvert() first creates a new matrix and then copies the data from
4306:    the first matrix.  A related routine is MatCopy(), which copies the matrix
4307:    entries of one matrix to another already existing matrix context.

4309:    Cannot be used to convert a sequential matrix to parallel or parallel to sequential,
4310:    the MPI communicator of the generated matrix is always the same as the communicator
4311:    of the input matrix.

4313:    Level: intermediate

4315: .seealso: MatCopy(), MatDuplicate()
4316: @*/
4317: PetscErrorCode MatConvert(Mat mat, MatType newtype,MatReuse reuse,Mat *M)
4318: {
4320:   PetscBool      sametype,issame,flg,issymmetric,ishermitian;
4321:   char           convname[256],mtype[256];
4322:   Mat            B;

4328:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
4329:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
4330:   MatCheckPreallocated(mat,1);

4332:   PetscOptionsGetString(((PetscObject)mat)->options,((PetscObject)mat)->prefix,"-matconvert_type",mtype,sizeof(mtype),&flg);
4333:   if (flg) newtype = mtype;

4335:   PetscObjectTypeCompare((PetscObject)mat,newtype,&sametype);
4336:   PetscStrcmp(newtype,"same",&issame);
4337:   if ((reuse == MAT_INPLACE_MATRIX) && (mat != *M)) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"MAT_INPLACE_MATRIX requires same input and output matrix");
4338:   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");

4340:   if ((reuse == MAT_INPLACE_MATRIX) && (issame || sametype)) {
4341:     PetscInfo3(mat,"Early return for inplace %s %d %d\n",((PetscObject)mat)->type_name,sametype,issame);
4342:     return(0);
4343:   }

4345:   /* Cache Mat options because some converter use MatHeaderReplace  */
4346:   issymmetric = mat->symmetric;
4347:   ishermitian = mat->hermitian;

4349:   if ((sametype || issame) && (reuse==MAT_INITIAL_MATRIX) && mat->ops->duplicate) {
4350:     PetscInfo3(mat,"Calling duplicate for initial matrix %s %d %d\n",((PetscObject)mat)->type_name,sametype,issame);
4351:     (*mat->ops->duplicate)(mat,MAT_COPY_VALUES,M);
4352:   } else {
4353:     PetscErrorCode (*conv)(Mat, MatType,MatReuse,Mat*)=NULL;
4354:     const char     *prefix[3] = {"seq","mpi",""};
4355:     PetscInt       i;
4356:     /*
4357:        Order of precedence:
4358:        0) See if newtype is a superclass of the current matrix.
4359:        1) See if a specialized converter is known to the current matrix.
4360:        2) See if a specialized converter is known to the desired matrix class.
4361:        3) See if a good general converter is registered for the desired class
4362:           (as of 6/27/03 only MATMPIADJ falls into this category).
4363:        4) See if a good general converter is known for the current matrix.
4364:        5) Use a really basic converter.
4365:     */

4367:     /* 0) See if newtype is a superclass of the current matrix.
4368:           i.e mat is mpiaij and newtype is aij */
4369:     for (i=0; i<2; i++) {
4370:       PetscStrncpy(convname,prefix[i],sizeof(convname));
4371:       PetscStrlcat(convname,newtype,sizeof(convname));
4372:       PetscStrcmp(convname,((PetscObject)mat)->type_name,&flg);
4373:       PetscInfo3(mat,"Check superclass %s %s -> %d\n",convname,((PetscObject)mat)->type_name,flg);
4374:       if (flg) {
4375:         if (reuse == MAT_INPLACE_MATRIX) {
4376:           PetscInfo(mat,"Early return\n");
4377:           return(0);
4378:         } else if (reuse == MAT_INITIAL_MATRIX && mat->ops->duplicate) {
4379:           PetscInfo(mat,"Calling MatDuplicate\n");
4380:           (*mat->ops->duplicate)(mat,MAT_COPY_VALUES,M);
4381:           return(0);
4382:         } else if (reuse == MAT_REUSE_MATRIX && mat->ops->copy) {
4383:           PetscInfo(mat,"Calling MatCopy\n");
4384:           MatCopy(mat,*M,SAME_NONZERO_PATTERN);
4385:           return(0);
4386:         }
4387:       }
4388:     }
4389:     /* 1) See if a specialized converter is known to the current matrix and the desired class */
4390:     for (i=0; i<3; i++) {
4391:       PetscStrncpy(convname,"MatConvert_",sizeof(convname));
4392:       PetscStrlcat(convname,((PetscObject)mat)->type_name,sizeof(convname));
4393:       PetscStrlcat(convname,"_",sizeof(convname));
4394:       PetscStrlcat(convname,prefix[i],sizeof(convname));
4395:       PetscStrlcat(convname,issame ? ((PetscObject)mat)->type_name : newtype,sizeof(convname));
4396:       PetscStrlcat(convname,"_C",sizeof(convname));
4397:       PetscObjectQueryFunction((PetscObject)mat,convname,&conv);
4398:       PetscInfo3(mat,"Check specialized (1) %s (%s) -> %d\n",convname,((PetscObject)mat)->type_name,!!conv);
4399:       if (conv) goto foundconv;
4400:     }

4402:     /* 2)  See if a specialized converter is known to the desired matrix class. */
4403:     MatCreate(PetscObjectComm((PetscObject)mat),&B);
4404:     MatSetSizes(B,mat->rmap->n,mat->cmap->n,mat->rmap->N,mat->cmap->N);
4405:     MatSetType(B,newtype);
4406:     for (i=0; i<3; i++) {
4407:       PetscStrncpy(convname,"MatConvert_",sizeof(convname));
4408:       PetscStrlcat(convname,((PetscObject)mat)->type_name,sizeof(convname));
4409:       PetscStrlcat(convname,"_",sizeof(convname));
4410:       PetscStrlcat(convname,prefix[i],sizeof(convname));
4411:       PetscStrlcat(convname,newtype,sizeof(convname));
4412:       PetscStrlcat(convname,"_C",sizeof(convname));
4413:       PetscObjectQueryFunction((PetscObject)B,convname,&conv);
4414:       PetscInfo3(mat,"Check specialized (2) %s (%s) -> %d\n",convname,((PetscObject)B)->type_name,!!conv);
4415:       if (conv) {
4416:         MatDestroy(&B);
4417:         goto foundconv;
4418:       }
4419:     }

4421:     /* 3) See if a good general converter is registered for the desired class */
4422:     conv = B->ops->convertfrom;
4423:     PetscInfo2(mat,"Check convertfrom (%s) -> %d\n",((PetscObject)B)->type_name,!!conv);
4424:     MatDestroy(&B);
4425:     if (conv) goto foundconv;

4427:     /* 4) See if a good general converter is known for the current matrix */
4428:     if (mat->ops->convert) conv = mat->ops->convert;

4430:     PetscInfo2(mat,"Check general convert (%s) -> %d\n",((PetscObject)mat)->type_name,!!conv);
4431:     if (conv) goto foundconv;

4433:     /* 5) Use a really basic converter. */
4434:     PetscInfo(mat,"Using MatConvert_Basic\n");
4435:     conv = MatConvert_Basic;

4437: foundconv:
4438:     PetscLogEventBegin(MAT_Convert,mat,0,0,0);
4439:     (*conv)(mat,newtype,reuse,M);
4440:     if (mat->rmap->mapping && mat->cmap->mapping && !(*M)->rmap->mapping && !(*M)->cmap->mapping) {
4441:       /* the block sizes must be same if the mappings are copied over */
4442:       (*M)->rmap->bs = mat->rmap->bs;
4443:       (*M)->cmap->bs = mat->cmap->bs;
4444:       PetscObjectReference((PetscObject)mat->rmap->mapping);
4445:       PetscObjectReference((PetscObject)mat->cmap->mapping);
4446:       (*M)->rmap->mapping = mat->rmap->mapping;
4447:       (*M)->cmap->mapping = mat->cmap->mapping;
4448:     }
4449:     (*M)->stencil.dim = mat->stencil.dim;
4450:     (*M)->stencil.noc = mat->stencil.noc;
4451:     for (i=0; i<=mat->stencil.dim; i++) {
4452:       (*M)->stencil.dims[i]   = mat->stencil.dims[i];
4453:       (*M)->stencil.starts[i] = mat->stencil.starts[i];
4454:     }
4455:     PetscLogEventEnd(MAT_Convert,mat,0,0,0);
4456:   }
4457:   PetscObjectStateIncrease((PetscObject)*M);

4459:   /* Copy Mat options */
4460:   if (issymmetric) {
4461:     MatSetOption(*M,MAT_SYMMETRIC,PETSC_TRUE);
4462:   }
4463:   if (ishermitian) {
4464:     MatSetOption(*M,MAT_HERMITIAN,PETSC_TRUE);
4465:   }
4466:   return(0);
4467: }

4469: /*@C
4470:    MatFactorGetSolverType - Returns name of the package providing the factorization routines

4472:    Not Collective

4474:    Input Parameter:
4475: .  mat - the matrix, must be a factored matrix

4477:    Output Parameter:
4478: .   type - the string name of the package (do not free this string)

4480:    Notes:
4481:       In Fortran you pass in a empty string and the package name will be copied into it.
4482:     (Make sure the string is long enough)

4484:    Level: intermediate

4486: .seealso: MatCopy(), MatDuplicate(), MatGetFactorAvailable(), MatGetFactor()
4487: @*/
4488: PetscErrorCode MatFactorGetSolverType(Mat mat, MatSolverType *type)
4489: {
4490:   PetscErrorCode ierr, (*conv)(Mat,MatSolverType*);

4495:   if (!mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
4496:   PetscObjectQueryFunction((PetscObject)mat,"MatFactorGetSolverType_C",&conv);
4497:   if (!conv) {
4498:     *type = MATSOLVERPETSC;
4499:   } else {
4500:     (*conv)(mat,type);
4501:   }
4502:   return(0);
4503: }

4505: typedef struct _MatSolverTypeForSpecifcType* MatSolverTypeForSpecifcType;
4506: struct _MatSolverTypeForSpecifcType {
4507:   MatType                        mtype;
4508:   PetscErrorCode                 (*createfactor[MAT_FACTOR_NUM_TYPES])(Mat,MatFactorType,Mat*);
4509:   MatSolverTypeForSpecifcType next;
4510: };

4512: typedef struct _MatSolverTypeHolder* MatSolverTypeHolder;
4513: struct _MatSolverTypeHolder {
4514:   char                        *name;
4515:   MatSolverTypeForSpecifcType handlers;
4516:   MatSolverTypeHolder         next;
4517: };

4519: static MatSolverTypeHolder MatSolverTypeHolders = NULL;

4521: /*@C
4522:    MatSolverTypeRegister - Registers a MatSolverType that works for a particular matrix type

4524:    Input Parameters:
4525: +    package - name of the package, for example petsc or superlu
4526: .    mtype - the matrix type that works with this package
4527: .    ftype - the type of factorization supported by the package
4528: -    createfactor - routine that will create the factored matrix ready to be used

4530:     Level: intermediate

4532: .seealso: MatCopy(), MatDuplicate(), MatGetFactorAvailable(), MatGetFactor()
4533: @*/
4534: PetscErrorCode MatSolverTypeRegister(MatSolverType package,MatType mtype,MatFactorType ftype,PetscErrorCode (*createfactor)(Mat,MatFactorType,Mat*))
4535: {
4536:   PetscErrorCode              ierr;
4537:   MatSolverTypeHolder         next = MatSolverTypeHolders,prev = NULL;
4538:   PetscBool                   flg;
4539:   MatSolverTypeForSpecifcType inext,iprev = NULL;

4542:   MatInitializePackage();
4543:   if (!next) {
4544:     PetscNew(&MatSolverTypeHolders);
4545:     PetscStrallocpy(package,&MatSolverTypeHolders->name);
4546:     PetscNew(&MatSolverTypeHolders->handlers);
4547:     PetscStrallocpy(mtype,(char **)&MatSolverTypeHolders->handlers->mtype);
4548:     MatSolverTypeHolders->handlers->createfactor[(int)ftype-1] = createfactor;
4549:     return(0);
4550:   }
4551:   while (next) {
4552:     PetscStrcasecmp(package,next->name,&flg);
4553:     if (flg) {
4554:       if (!next->handlers) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatSolverTypeHolder is missing handlers");
4555:       inext = next->handlers;
4556:       while (inext) {
4557:         PetscStrcasecmp(mtype,inext->mtype,&flg);
4558:         if (flg) {
4559:           inext->createfactor[(int)ftype-1] = createfactor;
4560:           return(0);
4561:         }
4562:         iprev = inext;
4563:         inext = inext->next;
4564:       }
4565:       PetscNew(&iprev->next);
4566:       PetscStrallocpy(mtype,(char **)&iprev->next->mtype);
4567:       iprev->next->createfactor[(int)ftype-1] = createfactor;
4568:       return(0);
4569:     }
4570:     prev = next;
4571:     next = next->next;
4572:   }
4573:   PetscNew(&prev->next);
4574:   PetscStrallocpy(package,&prev->next->name);
4575:   PetscNew(&prev->next->handlers);
4576:   PetscStrallocpy(mtype,(char **)&prev->next->handlers->mtype);
4577:   prev->next->handlers->createfactor[(int)ftype-1] = createfactor;
4578:   return(0);
4579: }

4581: /*@C
4582:    MatSolveTypeGet - Gets the function that creates the factor matrix if it exist

4584:    Input Parameters:
4585: +    type - name of the package, for example petsc or superlu
4586: .    ftype - the type of factorization supported by the type
4587: -    mtype - the matrix type that works with this type

4589:    Output Parameters:
4590: +   foundtype - PETSC_TRUE if the type was registered
4591: .   foundmtype - PETSC_TRUE if the type supports the requested mtype
4592: -   createfactor - routine that will create the factored matrix ready to be used or NULL if not found

4594:     Level: intermediate

4596: .seealso: MatCopy(), MatDuplicate(), MatGetFactorAvailable(), MatSolverTypeRegister(), MatGetFactor()
4597: @*/
4598: PetscErrorCode MatSolverTypeGet(MatSolverType type,MatType mtype,MatFactorType ftype,PetscBool *foundtype,PetscBool *foundmtype,PetscErrorCode (**createfactor)(Mat,MatFactorType,Mat*))
4599: {
4600:   PetscErrorCode              ierr;
4601:   MatSolverTypeHolder         next = MatSolverTypeHolders;
4602:   PetscBool                   flg;
4603:   MatSolverTypeForSpecifcType inext;

4606:   if (foundtype) *foundtype = PETSC_FALSE;
4607:   if (foundmtype)   *foundmtype   = PETSC_FALSE;
4608:   if (createfactor) *createfactor    = NULL;

4610:   if (type) {
4611:     while (next) {
4612:       PetscStrcasecmp(type,next->name,&flg);
4613:       if (flg) {
4614:         if (foundtype) *foundtype = PETSC_TRUE;
4615:         inext = next->handlers;
4616:         while (inext) {
4617:           PetscStrbeginswith(mtype,inext->mtype,&flg);
4618:           if (flg) {
4619:             if (foundmtype) *foundmtype = PETSC_TRUE;
4620:             if (createfactor)  *createfactor  = inext->createfactor[(int)ftype-1];
4621:             return(0);
4622:           }
4623:           inext = inext->next;
4624:         }
4625:       }
4626:       next = next->next;
4627:     }
4628:   } else {
4629:     while (next) {
4630:       inext = next->handlers;
4631:       while (inext) {
4632:         PetscStrcmp(mtype,inext->mtype,&flg);
4633:         if (flg && inext->createfactor[(int)ftype-1]) {
4634:           if (foundtype) *foundtype = PETSC_TRUE;
4635:           if (foundmtype)   *foundmtype   = PETSC_TRUE;
4636:           if (createfactor) *createfactor = inext->createfactor[(int)ftype-1];
4637:           return(0);
4638:         }
4639:         inext = inext->next;
4640:       }
4641:       next = next->next;
4642:     }
4643:     /* try with base classes inext->mtype */
4644:     next = MatSolverTypeHolders;
4645:     while (next) {
4646:       inext = next->handlers;
4647:       while (inext) {
4648:         PetscStrbeginswith(mtype,inext->mtype,&flg);
4649:         if (flg && inext->createfactor[(int)ftype-1]) {
4650:           if (foundtype) *foundtype = PETSC_TRUE;
4651:           if (foundmtype)   *foundmtype   = PETSC_TRUE;
4652:           if (createfactor) *createfactor = inext->createfactor[(int)ftype-1];
4653:           return(0);
4654:         }
4655:         inext = inext->next;
4656:       }
4657:       next = next->next;
4658:     }
4659:   }
4660:   return(0);
4661: }

4663: PetscErrorCode MatSolverTypeDestroy(void)
4664: {
4665:   PetscErrorCode              ierr;
4666:   MatSolverTypeHolder         next = MatSolverTypeHolders,prev;
4667:   MatSolverTypeForSpecifcType inext,iprev;

4670:   while (next) {
4671:     PetscFree(next->name);
4672:     inext = next->handlers;
4673:     while (inext) {
4674:       PetscFree(inext->mtype);
4675:       iprev = inext;
4676:       inext = inext->next;
4677:       PetscFree(iprev);
4678:     }
4679:     prev = next;
4680:     next = next->next;
4681:     PetscFree(prev);
4682:   }
4683:   MatSolverTypeHolders = NULL;
4684:   return(0);
4685: }

4687: /*@C
4688:    MatFactorGetUseOrdering - Indicates if the factorization uses the ordering provided in MatLUFactorSymbolic(), MatCholeskyFactorSymbolic()

4690:    Logically Collective on Mat

4692:    Input Parameters:
4693: .  mat - the matrix

4695:    Output Parameters:
4696: .  flg - PETSC_TRUE if uses the ordering

4698:    Notes:
4699:       Most internal PETSc factorizations use the ordering past to the factorization routine but external
4700:       packages do no, thus we want to skip the ordering when it is not needed.

4702:    Level: developer

4704: .seealso: MatCopy(), MatDuplicate(), MatGetFactorAvailable(), MatGetFactor(), MatLUFactorSymbolic(), MatCholeskyFactorSymbolic()
4705: @*/
4706: PetscErrorCode MatFactorGetUseOrdering(Mat mat, PetscBool *flg)
4707: {
4709:   *flg = mat->useordering;
4710:   return(0);
4711: }

4713: /*@C
4714:    MatGetFactor - Returns a matrix suitable to calls to MatXXFactorSymbolic()

4716:    Collective on Mat

4718:    Input Parameters:
4719: +  mat - the matrix
4720: .  type - name of solver type, for example, superlu, petsc (to use PETSc's default)
4721: -  ftype - factor type, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ICC, MAT_FACTOR_ILU,

4723:    Output Parameters:
4724: .  f - the factor matrix used with MatXXFactorSymbolic() calls

4726:    Notes:
4727:       Some PETSc matrix formats have alternative solvers available that are contained in alternative packages
4728:      such as pastix, superlu, mumps etc.

4730:       PETSc must have been ./configure to use the external solver, using the option --download-package

4732:    Developer Notes:
4733:       This should actually be called MatCreateFactor() since it creates a new factor object

4735:    Level: intermediate

4737: .seealso: MatCopy(), MatDuplicate(), MatGetFactorAvailable(), MatFactorGetUseOrdering(), MatSolverTypeRegister()
4738: @*/
4739: PetscErrorCode MatGetFactor(Mat mat, MatSolverType type,MatFactorType ftype,Mat *f)
4740: {
4741:   PetscErrorCode ierr,(*conv)(Mat,MatFactorType,Mat*);
4742:   PetscBool      foundtype,foundmtype;


4748:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
4749:   MatCheckPreallocated(mat,1);

4751:   MatSolverTypeGet(type,((PetscObject)mat)->type_name,ftype,&foundtype,&foundmtype,&conv);
4752:   if (!foundtype) {
4753:     if (type) {
4754:       SETERRQ4(PetscObjectComm((PetscObject)mat),PETSC_ERR_MISSING_FACTOR,"Could not locate solver type %s for factorization type %s and matrix type %s. Perhaps you must ./configure with --download-%s",type,MatFactorTypes[ftype],((PetscObject)mat)->type_name,type);
4755:     } else {
4756:       SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_MISSING_FACTOR,"Could not locate a solver type for factorization type %s and matrix type %s.",MatFactorTypes[ftype],((PetscObject)mat)->type_name);
4757:     }
4758:   }
4759:   if (!foundmtype) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_MISSING_FACTOR,"MatSolverType %s does not support matrix type %s",type,((PetscObject)mat)->type_name);
4760:   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);

4762:   (*conv)(mat,ftype,f);
4763:   return(0);
4764: }

4766: /*@C
4767:    MatGetFactorAvailable - Returns a a flag if matrix supports particular type and factor type

4769:    Not Collective

4771:    Input Parameters:
4772: +  mat - the matrix
4773: .  type - name of solver type, for example, superlu, petsc (to use PETSc's default)
4774: -  ftype - factor type, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ICC, MAT_FACTOR_ILU,

4776:    Output Parameter:
4777: .    flg - PETSC_TRUE if the factorization is available

4779:    Notes:
4780:       Some PETSc matrix formats have alternative solvers available that are contained in alternative packages
4781:      such as pastix, superlu, mumps etc.

4783:       PETSc must have been ./configure to use the external solver, using the option --download-package

4785:    Developer Notes:
4786:       This should actually be called MatCreateFactorAvailable() since MatGetFactor() creates a new factor object

4788:    Level: intermediate

4790: .seealso: MatCopy(), MatDuplicate(), MatGetFactor(), MatSolverTypeRegister()
4791: @*/
4792: PetscErrorCode MatGetFactorAvailable(Mat mat, MatSolverType type,MatFactorType ftype,PetscBool  *flg)
4793: {
4794:   PetscErrorCode ierr, (*gconv)(Mat,MatFactorType,Mat*);


4800:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
4801:   MatCheckPreallocated(mat,1);

4803:   *flg = PETSC_FALSE;
4804:   MatSolverTypeGet(type,((PetscObject)mat)->type_name,ftype,NULL,NULL,&gconv);
4805:   if (gconv) {
4806:     *flg = PETSC_TRUE;
4807:   }
4808:   return(0);
4809: }

4811: #include <petscdmtypes.h>

4813: /*@
4814:    MatDuplicate - Duplicates a matrix including the non-zero structure.

4816:    Collective on Mat

4818:    Input Parameters:
4819: +  mat - the matrix
4820: -  op - One of MAT_DO_NOT_COPY_VALUES, MAT_COPY_VALUES, or MAT_SHARE_NONZERO_PATTERN.
4821:         See the manual page for MatDuplicateOption for an explanation of these options.

4823:    Output Parameter:
4824: .  M - pointer to place new matrix

4826:    Level: intermediate

4828:    Notes:
4829:     You cannot change the nonzero pattern for the parent or child matrix if you use MAT_SHARE_NONZERO_PATTERN.
4830:     When original mat is a product of matrix operation, e.g., an output of MatMatMult() or MatCreateSubMatrix(), only the simple matrix data structure of mat is duplicated and the internal data structures created for the reuse of previous matrix operations are not duplicated. User should not use MatDuplicate() to create new matrix M if M is intended to be reused as the product of matrix operation.

4832: .seealso: MatCopy(), MatConvert(), MatDuplicateOption
4833: @*/
4834: PetscErrorCode MatDuplicate(Mat mat,MatDuplicateOption op,Mat *M)
4835: {
4837:   Mat            B;
4838:   PetscInt       i;
4839:   DM             dm;
4840:   void           (*viewf)(void);

4846:   if (op == MAT_COPY_VALUES && !mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"MAT_COPY_VALUES not allowed for unassembled matrix");
4847:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
4848:   MatCheckPreallocated(mat,1);

4850:   *M = NULL;
4851:   if (!mat->ops->duplicate) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Not written for matrix type %s\n",((PetscObject)mat)->type_name);
4852:   PetscLogEventBegin(MAT_Convert,mat,0,0,0);
4853:   (*mat->ops->duplicate)(mat,op,M);
4854:   B    = *M;

4856:   MatGetOperation(mat,MATOP_VIEW,&viewf);
4857:   if (viewf) {
4858:     MatSetOperation(B,MATOP_VIEW,viewf);
4859:   }

4861:   B->stencil.dim = mat->stencil.dim;
4862:   B->stencil.noc = mat->stencil.noc;
4863:   for (i=0; i<=mat->stencil.dim; i++) {
4864:     B->stencil.dims[i]   = mat->stencil.dims[i];
4865:     B->stencil.starts[i] = mat->stencil.starts[i];
4866:   }

4868:   B->nooffproczerorows = mat->nooffproczerorows;
4869:   B->nooffprocentries  = mat->nooffprocentries;

4871:   PetscObjectQuery((PetscObject) mat, "__PETSc_dm", (PetscObject*) &dm);
4872:   if (dm) {
4873:     PetscObjectCompose((PetscObject) B, "__PETSc_dm", (PetscObject) dm);
4874:   }
4875:   PetscLogEventEnd(MAT_Convert,mat,0,0,0);
4876:   PetscObjectStateIncrease((PetscObject)B);
4877:   return(0);
4878: }

4880: /*@
4881:    MatGetDiagonal - Gets the diagonal of a matrix.

4883:    Logically Collective on Mat

4885:    Input Parameters:
4886: +  mat - the matrix
4887: -  v - the vector for storing the diagonal

4889:    Output Parameter:
4890: .  v - the diagonal of the matrix

4892:    Level: intermediate

4894:    Note:
4895:    Currently only correct in parallel for square matrices.

4897: .seealso: MatGetRow(), MatCreateSubMatrices(), MatCreateSubMatrix(), MatGetRowMaxAbs()
4898: @*/
4899: PetscErrorCode MatGetDiagonal(Mat mat,Vec v)
4900: {

4907:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
4908:   if (!mat->ops->getdiagonal) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
4909:   MatCheckPreallocated(mat,1);

4911:   (*mat->ops->getdiagonal)(mat,v);
4912:   PetscObjectStateIncrease((PetscObject)v);
4913:   return(0);
4914: }

4916: /*@C
4917:    MatGetRowMin - Gets the minimum value (of the real part) of each
4918:         row of the matrix

4920:    Logically Collective on Mat

4922:    Input Parameters:
4923: .  mat - the matrix

4925:    Output Parameter:
4926: +  v - the vector for storing the maximums
4927: -  idx - the indices of the column found for each row (optional)

4929:    Level: intermediate

4931:    Notes:
4932:     The result of this call are the same as if one converted the matrix to dense format
4933:       and found the minimum value in each row (i.e. the implicit zeros are counted as zeros).

4935:     This code is only implemented for a couple of matrix formats.

4937: .seealso: MatGetDiagonal(), MatCreateSubMatrices(), MatCreateSubMatrix(), MatGetRowMaxAbs(),
4938:           MatGetRowMax()
4939: @*/
4940: PetscErrorCode MatGetRowMin(Mat mat,Vec v,PetscInt idx[])
4941: {

4948:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");

4950:   if (!mat->cmap->N) {
4951:     VecSet(v,PETSC_MAX_REAL);
4952:     if (idx) {
4953:       PetscInt i,m = mat->rmap->n;
4954:       for (i=0; i<m; i++) idx[i] = -1;
4955:     }
4956:   } else {
4957:     if (!mat->ops->getrowmin) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
4958:     MatCheckPreallocated(mat,1);
4959:   }
4960:   (*mat->ops->getrowmin)(mat,v,idx);
4961:   PetscObjectStateIncrease((PetscObject)v);
4962:   return(0);
4963: }

4965: /*@C
4966:    MatGetRowMinAbs - Gets the minimum value (in absolute value) of each
4967:         row of the matrix

4969:    Logically Collective on Mat

4971:    Input Parameters:
4972: .  mat - the matrix

4974:    Output Parameter:
4975: +  v - the vector for storing the minimums
4976: -  idx - the indices of the column found for each row (or NULL if not needed)

4978:    Level: intermediate

4980:    Notes:
4981:     if a row is completely empty or has only 0.0 values then the idx[] value for that
4982:     row is 0 (the first column).

4984:     This code is only implemented for a couple of matrix formats.

4986: .seealso: MatGetDiagonal(), MatCreateSubMatrices(), MatCreateSubMatrix(), MatGetRowMax(), MatGetRowMaxAbs(), MatGetRowMin()
4987: @*/
4988: PetscErrorCode MatGetRowMinAbs(Mat mat,Vec v,PetscInt idx[])
4989: {

4996:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
4997:   if (mat->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");

4999:   if (!mat->cmap->N) {
5000:     VecSet(v,0.0);
5001:     if (idx) {
5002:       PetscInt i,m = mat->rmap->n;
5003:       for (i=0; i<m; i++) idx[i] = -1;
5004:     }
5005:   } else {
5006:     if (!mat->ops->getrowminabs) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
5007:     MatCheckPreallocated(mat,1);
5008:     if (idx) {PetscArrayzero(idx,mat->rmap->n);}
5009:     (*mat->ops->getrowminabs)(mat,v,idx);
5010:   }
5011:   PetscObjectStateIncrease((PetscObject)v);
5012:   return(0);
5013: }

5015: /*@C
5016:    MatGetRowMax - Gets the maximum value (of the real part) of each
5017:         row of the matrix

5019:    Logically Collective on Mat

5021:    Input Parameters:
5022: .  mat - the matrix

5024:    Output Parameter:
5025: +  v - the vector for storing the maximums
5026: -  idx - the indices of the column found for each row (optional)

5028:    Level: intermediate

5030:    Notes:
5031:     The result of this call are the same as if one converted the matrix to dense format
5032:       and found the minimum value in each row (i.e. the implicit zeros are counted as zeros).

5034:     This code is only implemented for a couple of matrix formats.

5036: .seealso: MatGetDiagonal(), MatCreateSubMatrices(), MatCreateSubMatrix(), MatGetRowMaxAbs(), MatGetRowMin()
5037: @*/
5038: PetscErrorCode MatGetRowMax(Mat mat,Vec v,PetscInt idx[])
5039: {

5046:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");

5048:   if (!mat->cmap->N) {
5049:     VecSet(v,PETSC_MIN_REAL);
5050:     if (idx) {
5051:       PetscInt i,m = mat->rmap->n;
5052:       for (i=0; i<m; i++) idx[i] = -1;
5053:     }
5054:   } else {
5055:     if (!mat->ops->getrowmax) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
5056:     MatCheckPreallocated(mat,1);
5057:     (*mat->ops->getrowmax)(mat,v,idx);
5058:   }
5059:   PetscObjectStateIncrease((PetscObject)v);
5060:   return(0);
5061: }

5063: /*@C
5064:    MatGetRowMaxAbs - Gets the maximum value (in absolute value) of each
5065:         row of the matrix

5067:    Logically Collective on Mat

5069:    Input Parameters:
5070: .  mat - the matrix

5072:    Output Parameter:
5073: +  v - the vector for storing the maximums
5074: -  idx - the indices of the column found for each row (or NULL if not needed)

5076:    Level: intermediate

5078:    Notes:
5079:     if a row is completely empty or has only 0.0 values then the idx[] value for that
5080:     row is 0 (the first column).

5082:     This code is only implemented for a couple of matrix formats.

5084: .seealso: MatGetDiagonal(), MatCreateSubMatrices(), MatCreateSubMatrix(), MatGetRowMax(), MatGetRowMin()
5085: @*/
5086: PetscErrorCode MatGetRowMaxAbs(Mat mat,Vec v,PetscInt idx[])
5087: {

5094:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");

5096:   if (!mat->cmap->N) {
5097:     VecSet(v,0.0);
5098:     if (idx) {
5099:       PetscInt i,m = mat->rmap->n;
5100:       for (i=0; i<m; i++) idx[i] = -1;
5101:     }
5102:   } else {
5103:     if (!mat->ops->getrowmaxabs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
5104:     MatCheckPreallocated(mat,1);
5105:     if (idx) {PetscArrayzero(idx,mat->rmap->n);}
5106:     (*mat->ops->getrowmaxabs)(mat,v,idx);
5107:   }
5108:   PetscObjectStateIncrease((PetscObject)v);
5109:   return(0);
5110: }

5112: /*@
5113:    MatGetRowSum - Gets the sum of each row of the matrix

5115:    Logically or Neighborhood Collective on Mat

5117:    Input Parameters:
5118: .  mat - the matrix

5120:    Output Parameter:
5121: .  v - the vector for storing the sum of rows

5123:    Level: intermediate

5125:    Notes:
5126:     This code is slow since it is not currently specialized for different formats

5128: .seealso: MatGetDiagonal(), MatCreateSubMatrices(), MatCreateSubMatrix(), MatGetRowMax(), MatGetRowMin()
5129: @*/
5130: PetscErrorCode MatGetRowSum(Mat mat, Vec v)
5131: {
5132:   Vec            ones;

5139:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
5140:   MatCheckPreallocated(mat,1);
5141:   MatCreateVecs(mat,&ones,NULL);
5142:   VecSet(ones,1.);
5143:   MatMult(mat,ones,v);
5144:   VecDestroy(&ones);
5145:   return(0);
5146: }

5148: /*@
5149:    MatTranspose - Computes an in-place or out-of-place transpose of a matrix.

5151:    Collective on Mat

5153:    Input Parameters:
5154: +  mat - the matrix to transpose
5155: -  reuse - either MAT_INITIAL_MATRIX, MAT_REUSE_MATRIX, or MAT_INPLACE_MATRIX

5157:    Output Parameter:
5158: .  B - the transpose

5160:    Notes:
5161:      If you use MAT_INPLACE_MATRIX then you must pass in &mat for B

5163:      MAT_REUSE_MATRIX causes the B matrix from a previous call to this function with MAT_INITIAL_MATRIX to be used

5165:      Consider using MatCreateTranspose() instead if you only need a matrix that behaves like the transpose, but don't need the storage to be changed.

5167:    Level: intermediate

5169: .seealso: MatMultTranspose(), MatMultTransposeAdd(), MatIsTranspose(), MatReuse
5170: @*/
5171: PetscErrorCode MatTranspose(Mat mat,MatReuse reuse,Mat *B)
5172: {

5178:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
5179:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
5180:   if (!mat->ops->transpose) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
5181:   if (reuse == MAT_INPLACE_MATRIX && mat != *B) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"MAT_INPLACE_MATRIX requires last matrix to match first");
5182:   if (reuse == MAT_REUSE_MATRIX && mat == *B) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Perhaps you mean MAT_INPLACE_MATRIX");
5183:   MatCheckPreallocated(mat,1);

5185:   PetscLogEventBegin(MAT_Transpose,mat,0,0,0);
5186:   (*mat->ops->transpose)(mat,reuse,B);
5187:   PetscLogEventEnd(MAT_Transpose,mat,0,0,0);
5188:   if (B) {PetscObjectStateIncrease((PetscObject)*B);}
5189:   return(0);
5190: }

5192: /*@
5193:    MatIsTranspose - Test whether a matrix is another one's transpose,
5194:         or its own, in which case it tests symmetry.

5196:    Collective on Mat

5198:    Input Parameter:
5199: +  A - the matrix to test
5200: -  B - the matrix to test against, this can equal the first parameter

5202:    Output Parameters:
5203: .  flg - the result

5205:    Notes:
5206:    Only available for SeqAIJ/MPIAIJ matrices. The sequential algorithm
5207:    has a running time of the order of the number of nonzeros; the parallel
5208:    test involves parallel copies of the block-offdiagonal parts of the matrix.

5210:    Level: intermediate

5212: .seealso: MatTranspose(), MatIsSymmetric(), MatIsHermitian()
5213: @*/
5214: PetscErrorCode MatIsTranspose(Mat A,Mat B,PetscReal tol,PetscBool  *flg)
5215: {
5216:   PetscErrorCode ierr,(*f)(Mat,Mat,PetscReal,PetscBool*),(*g)(Mat,Mat,PetscReal,PetscBool*);

5222:   PetscObjectQueryFunction((PetscObject)A,"MatIsTranspose_C",&f);
5223:   PetscObjectQueryFunction((PetscObject)B,"MatIsTranspose_C",&g);
5224:   *flg = PETSC_FALSE;
5225:   if (f && g) {
5226:     if (f == g) {
5227:       (*f)(A,B,tol,flg);
5228:     } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_NOTSAMETYPE,"Matrices do not have the same comparator for symmetry test");
5229:   } else {
5230:     MatType mattype;
5231:     if (!f) {
5232:       MatGetType(A,&mattype);
5233:     } else {
5234:       MatGetType(B,&mattype);
5235:     }
5236:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix of type %s does not support checking for transpose",mattype);
5237:   }
5238:   return(0);
5239: }

5241: /*@
5242:    MatHermitianTranspose - Computes an in-place or out-of-place transpose of a matrix in complex conjugate.

5244:    Collective on Mat

5246:    Input Parameters:
5247: +  mat - the matrix to transpose and complex conjugate
5248: -  reuse - either MAT_INITIAL_MATRIX, MAT_REUSE_MATRIX, or MAT_INPLACE_MATRIX

5250:    Output Parameter:
5251: .  B - the Hermitian

5253:    Level: intermediate

5255: .seealso: MatTranspose(), MatMultTranspose(), MatMultTransposeAdd(), MatIsTranspose(), MatReuse
5256: @*/
5257: PetscErrorCode MatHermitianTranspose(Mat mat,MatReuse reuse,Mat *B)
5258: {

5262:   MatTranspose(mat,reuse,B);
5263: #if defined(PETSC_USE_COMPLEX)
5264:   MatConjugate(*B);
5265: #endif
5266:   return(0);
5267: }

5269: /*@
5270:    MatIsHermitianTranspose - Test whether a matrix is another one's Hermitian transpose,

5272:    Collective on Mat

5274:    Input Parameter:
5275: +  A - the matrix to test
5276: -  B - the matrix to test against, this can equal the first parameter

5278:    Output Parameters:
5279: .  flg - the result

5281:    Notes:
5282:    Only available for SeqAIJ/MPIAIJ matrices. The sequential algorithm
5283:    has a running time of the order of the number of nonzeros; the parallel
5284:    test involves parallel copies of the block-offdiagonal parts of the matrix.

5286:    Level: intermediate

5288: .seealso: MatTranspose(), MatIsSymmetric(), MatIsHermitian(), MatIsTranspose()
5289: @*/
5290: PetscErrorCode MatIsHermitianTranspose(Mat A,Mat B,PetscReal tol,PetscBool  *flg)
5291: {
5292:   PetscErrorCode ierr,(*f)(Mat,Mat,PetscReal,PetscBool*),(*g)(Mat,Mat,PetscReal,PetscBool*);

5298:   PetscObjectQueryFunction((PetscObject)A,"MatIsHermitianTranspose_C",&f);
5299:   PetscObjectQueryFunction((PetscObject)B,"MatIsHermitianTranspose_C",&g);
5300:   if (f && g) {
5301:     if (f==g) {
5302:       (*f)(A,B,tol,flg);
5303:     } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_NOTSAMETYPE,"Matrices do not have the same comparator for Hermitian test");
5304:   }
5305:   return(0);
5306: }

5308: /*@
5309:    MatPermute - Creates a new matrix with rows and columns permuted from the
5310:    original.

5312:    Collective on Mat

5314:    Input Parameters:
5315: +  mat - the matrix to permute
5316: .  row - row permutation, each processor supplies only the permutation for its rows
5317: -  col - column permutation, each processor supplies only the permutation for its columns

5319:    Output Parameters:
5320: .  B - the permuted matrix

5322:    Level: advanced

5324:    Note:
5325:    The index sets map from row/col of permuted matrix to row/col of original matrix.
5326:    The index sets should be on the same communicator as Mat and have the same local sizes.

5328:    Developer Note:
5329:      If you want to implement MatPermute for a matrix type, and your approach doesn't
5330:      exploit the fact that row and col are permutations, consider implementing the
5331:      more general MatCreateSubMatrix() instead.

5333: .seealso: MatGetOrdering(), ISAllGather()

5335: @*/
5336: PetscErrorCode MatPermute(Mat mat,IS row,IS col,Mat *B)
5337: {

5348:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
5349:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
5350:   if (!mat->ops->permute && !mat->ops->createsubmatrix) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatPermute not available for Mat type %s",((PetscObject)mat)->type_name);
5351:   MatCheckPreallocated(mat,1);

5353:   if (mat->ops->permute) {
5354:     (*mat->ops->permute)(mat,row,col,B);
5355:     PetscObjectStateIncrease((PetscObject)*B);
5356:   } else {
5357:     MatCreateSubMatrix(mat, row, col, MAT_INITIAL_MATRIX, B);
5358:   }
5359:   return(0);
5360: }

5362: /*@
5363:    MatEqual - Compares two matrices.

5365:    Collective on Mat

5367:    Input Parameters:
5368: +  A - the first matrix
5369: -  B - the second matrix

5371:    Output Parameter:
5372: .  flg - PETSC_TRUE if the matrices are equal; PETSC_FALSE otherwise.

5374:    Level: intermediate

5376: @*/
5377: PetscErrorCode MatEqual(Mat A,Mat B,PetscBool  *flg)
5378: {

5388:   MatCheckPreallocated(B,2);
5389:   if (!A->assembled) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
5390:   if (!B->assembled) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
5391:   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);
5392:   if (!A->ops->equal) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Mat type %s",((PetscObject)A)->type_name);
5393:   if (!B->ops->equal) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Mat type %s",((PetscObject)B)->type_name);
5394:   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);
5395:   MatCheckPreallocated(A,1);

5397:   (*A->ops->equal)(A,B,flg);
5398:   return(0);
5399: }

5401: /*@
5402:    MatDiagonalScale - Scales a matrix on the left and right by diagonal
5403:    matrices that are stored as vectors.  Either of the two scaling
5404:    matrices can be NULL.

5406:    Collective on Mat

5408:    Input Parameters:
5409: +  mat - the matrix to be scaled
5410: .  l - the left scaling vector (or NULL)
5411: -  r - the right scaling vector (or NULL)

5413:    Notes:
5414:    MatDiagonalScale() computes A = LAR, where
5415:    L = a diagonal matrix (stored as a vector), R = a diagonal matrix (stored as a vector)
5416:    The L scales the rows of the matrix, the R scales the columns of the matrix.

5418:    Level: intermediate


5421: .seealso: MatScale(), MatShift(), MatDiagonalSet()
5422: @*/
5423: PetscErrorCode MatDiagonalScale(Mat mat,Vec l,Vec r)
5424: {

5432:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
5433:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
5434:   MatCheckPreallocated(mat,1);
5435:   if (!l && !r) return(0);

5437:   if (!mat->ops->diagonalscale) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
5438:   PetscLogEventBegin(MAT_Scale,mat,0,0,0);
5439:   (*mat->ops->diagonalscale)(mat,l,r);
5440:   PetscLogEventEnd(MAT_Scale,mat,0,0,0);
5441:   PetscObjectStateIncrease((PetscObject)mat);
5442:   return(0);
5443: }

5445: /*@
5446:     MatScale - Scales all elements of a matrix by a given number.

5448:     Logically Collective on Mat

5450:     Input Parameters:
5451: +   mat - the matrix to be scaled
5452: -   a  - the scaling value

5454:     Output Parameter:
5455: .   mat - the scaled matrix

5457:     Level: intermediate

5459: .seealso: MatDiagonalScale()
5460: @*/
5461: PetscErrorCode MatScale(Mat mat,PetscScalar a)
5462: {

5468:   if (a != (PetscScalar)1.0 && !mat->ops->scale) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
5469:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
5470:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
5472:   MatCheckPreallocated(mat,1);

5474:   PetscLogEventBegin(MAT_Scale,mat,0,0,0);
5475:   if (a != (PetscScalar)1.0) {
5476:     (*mat->ops->scale)(mat,a);
5477:     PetscObjectStateIncrease((PetscObject)mat);
5478:   }
5479:   PetscLogEventEnd(MAT_Scale,mat,0,0,0);
5480:   return(0);
5481: }

5483: /*@
5484:    MatNorm - Calculates various norms of a matrix.

5486:    Collective on Mat

5488:    Input Parameters:
5489: +  mat - the matrix
5490: -  type - the type of norm, NORM_1, NORM_FROBENIUS, NORM_INFINITY

5492:    Output Parameters:
5493: .  nrm - the resulting norm

5495:    Level: intermediate

5497: @*/
5498: PetscErrorCode MatNorm(Mat mat,NormType type,PetscReal *nrm)
5499: {


5507:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
5508:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
5509:   if (!mat->ops->norm) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
5510:   MatCheckPreallocated(mat,1);

5512:   (*mat->ops->norm)(mat,type,nrm);
5513:   return(0);
5514: }

5516: /*
5517:      This variable is used to prevent counting of MatAssemblyBegin() that
5518:    are called from within a MatAssemblyEnd().
5519: */
5520: static PetscInt MatAssemblyEnd_InUse = 0;
5521: /*@
5522:    MatAssemblyBegin - Begins assembling the matrix.  This routine should
5523:    be called after completing all calls to MatSetValues().

5525:    Collective on Mat

5527:    Input Parameters:
5528: +  mat - the matrix
5529: -  type - type of assembly, either MAT_FLUSH_ASSEMBLY or MAT_FINAL_ASSEMBLY

5531:    Notes:
5532:    MatSetValues() generally caches the values.  The matrix is ready to
5533:    use only after MatAssemblyBegin() and MatAssemblyEnd() have been called.
5534:    Use MAT_FLUSH_ASSEMBLY when switching between ADD_VALUES and INSERT_VALUES
5535:    in MatSetValues(); use MAT_FINAL_ASSEMBLY for the final assembly before
5536:    using the matrix.

5538:    ALL processes that share a matrix MUST call MatAssemblyBegin() and MatAssemblyEnd() the SAME NUMBER of times, and each time with the
5539:    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
5540:    a global collective operation requring all processes that share the matrix.

5542:    Space for preallocated nonzeros that is not filled by a call to MatSetValues() or a related routine are compressed
5543:    out by assembly. If you intend to use that extra space on a subsequent assembly, be sure to insert explicit zeros
5544:    before MAT_FINAL_ASSEMBLY so the space is not compressed out.

5546:    Level: beginner

5548: .seealso: MatAssemblyEnd(), MatSetValues(), MatAssembled()
5549: @*/
5550: PetscErrorCode MatAssemblyBegin(Mat mat,MatAssemblyType type)
5551: {

5557:   MatCheckPreallocated(mat,1);
5558:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix.\nDid you forget to call MatSetUnfactored()?");
5559:   if (mat->assembled) {
5560:     mat->was_assembled = PETSC_TRUE;
5561:     mat->assembled     = PETSC_FALSE;
5562:   }

5564:   if (!MatAssemblyEnd_InUse) {
5565:     PetscLogEventBegin(MAT_AssemblyBegin,mat,0,0,0);
5566:     if (mat->ops->assemblybegin) {(*mat->ops->assemblybegin)(mat,type);}
5567:     PetscLogEventEnd(MAT_AssemblyBegin,mat,0,0,0);
5568:   } else if (mat->ops->assemblybegin) {
5569:     (*mat->ops->assemblybegin)(mat,type);
5570:   }
5571:   return(0);
5572: }

5574: /*@
5575:    MatAssembled - Indicates if a matrix has been assembled and is ready for
5576:      use; for example, in matrix-vector product.

5578:    Not Collective

5580:    Input Parameter:
5581: .  mat - the matrix

5583:    Output Parameter:
5584: .  assembled - PETSC_TRUE or PETSC_FALSE

5586:    Level: advanced

5588: .seealso: MatAssemblyEnd(), MatSetValues(), MatAssemblyBegin()
5589: @*/
5590: PetscErrorCode MatAssembled(Mat mat,PetscBool  *assembled)
5591: {
5595:   *assembled = mat->assembled;
5596:   return(0);
5597: }

5599: /*@
5600:    MatAssemblyEnd - Completes assembling the matrix.  This routine should
5601:    be called after MatAssemblyBegin().

5603:    Collective on Mat

5605:    Input Parameters:
5606: +  mat - the matrix
5607: -  type - type of assembly, either MAT_FLUSH_ASSEMBLY or MAT_FINAL_ASSEMBLY

5609:    Options Database Keys:
5610: +  -mat_view ::ascii_info - Prints info on matrix at conclusion of MatEndAssembly()
5611: .  -mat_view ::ascii_info_detail - Prints more detailed info
5612: .  -mat_view - Prints matrix in ASCII format
5613: .  -mat_view ::ascii_matlab - Prints matrix in Matlab format
5614: .  -mat_view draw - PetscDraws nonzero structure of matrix, using MatView() and PetscDrawOpenX().
5615: .  -display <name> - Sets display name (default is host)
5616: .  -draw_pause <sec> - Sets number of seconds to pause after display
5617: .  -mat_view socket - Sends matrix to socket, can be accessed from Matlab (See Users-Manual: ch_matlab)
5618: .  -viewer_socket_machine <machine> - Machine to use for socket
5619: .  -viewer_socket_port <port> - Port number to use for socket
5620: -  -mat_view binary:filename[:append] - Save matrix to file in binary format

5622:    Notes:
5623:    MatSetValues() generally caches the values.  The matrix is ready to
5624:    use only after MatAssemblyBegin() and MatAssemblyEnd() have been called.
5625:    Use MAT_FLUSH_ASSEMBLY when switching between ADD_VALUES and INSERT_VALUES
5626:    in MatSetValues(); use MAT_FINAL_ASSEMBLY for the final assembly before
5627:    using the matrix.

5629:    Space for preallocated nonzeros that is not filled by a call to MatSetValues() or a related routine are compressed
5630:    out by assembly. If you intend to use that extra space on a subsequent assembly, be sure to insert explicit zeros
5631:    before MAT_FINAL_ASSEMBLY so the space is not compressed out.

5633:    Level: beginner

5635: .seealso: MatAssemblyBegin(), MatSetValues(), PetscDrawOpenX(), PetscDrawCreate(), MatView(), MatAssembled(), PetscViewerSocketOpen()
5636: @*/
5637: PetscErrorCode MatAssemblyEnd(Mat mat,MatAssemblyType type)
5638: {
5639:   PetscErrorCode  ierr;
5640:   static PetscInt inassm = 0;
5641:   PetscBool       flg    = PETSC_FALSE;


5647:   inassm++;
5648:   MatAssemblyEnd_InUse++;
5649:   if (MatAssemblyEnd_InUse == 1) { /* Do the logging only the first time through */
5650:     PetscLogEventBegin(MAT_AssemblyEnd,mat,0,0,0);
5651:     if (mat->ops->assemblyend) {
5652:       (*mat->ops->assemblyend)(mat,type);
5653:     }
5654:     PetscLogEventEnd(MAT_AssemblyEnd,mat,0,0,0);
5655:   } else if (mat->ops->assemblyend) {
5656:     (*mat->ops->assemblyend)(mat,type);
5657:   }

5659:   /* Flush assembly is not a true assembly */
5660:   if (type != MAT_FLUSH_ASSEMBLY) {
5661:     mat->num_ass++;
5662:     mat->assembled        = PETSC_TRUE;
5663:     mat->ass_nonzerostate = mat->nonzerostate;
5664:   }

5666:   mat->insertmode = NOT_SET_VALUES;
5667:   MatAssemblyEnd_InUse--;
5668:   PetscObjectStateIncrease((PetscObject)mat);
5669:   if (!mat->symmetric_eternal) {
5670:     mat->symmetric_set              = PETSC_FALSE;
5671:     mat->hermitian_set              = PETSC_FALSE;
5672:     mat->structurally_symmetric_set = PETSC_FALSE;
5673:   }
5674:   if (inassm == 1 && type != MAT_FLUSH_ASSEMBLY) {
5675:     MatViewFromOptions(mat,NULL,"-mat_view");

5677:     if (mat->checksymmetryonassembly) {
5678:       MatIsSymmetric(mat,mat->checksymmetrytol,&flg);
5679:       if (flg) {
5680:         PetscPrintf(PetscObjectComm((PetscObject)mat),"Matrix is symmetric (tolerance %g)\n",(double)mat->checksymmetrytol);
5681:       } else {
5682:         PetscPrintf(PetscObjectComm((PetscObject)mat),"Matrix is not symmetric (tolerance %g)\n",(double)mat->checksymmetrytol);
5683:       }
5684:     }
5685:     if (mat->nullsp && mat->checknullspaceonassembly) {
5686:       MatNullSpaceTest(mat->nullsp,mat,NULL);
5687:     }
5688:   }
5689:   inassm--;
5690:   return(0);
5691: }

5693: /*@
5694:    MatSetOption - Sets a parameter option for a matrix. Some options
5695:    may be specific to certain storage formats.  Some options
5696:    determine how values will be inserted (or added). Sorted,
5697:    row-oriented input will generally assemble the fastest. The default
5698:    is row-oriented.

5700:    Logically Collective on Mat for certain operations, such as MAT_SPD, not collective for MAT_ROW_ORIENTED, see MatOption

5702:    Input Parameters:
5703: +  mat - the matrix
5704: .  option - the option, one of those listed below (and possibly others),
5705: -  flg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE)

5707:   Options Describing Matrix Structure:
5708: +    MAT_SPD - symmetric positive definite
5709: .    MAT_SYMMETRIC - symmetric in terms of both structure and value
5710: .    MAT_HERMITIAN - transpose is the complex conjugation
5711: .    MAT_STRUCTURALLY_SYMMETRIC - symmetric nonzero structure
5712: -    MAT_SYMMETRY_ETERNAL - if you would like the symmetry/Hermitian flag
5713:                             you set to be kept with all future use of the matrix
5714:                             including after MatAssemblyBegin/End() which could
5715:                             potentially change the symmetry structure, i.e. you
5716:                             KNOW the matrix will ALWAYS have the property you set.
5717:                             Note that setting this flag alone implies nothing about whether the matrix is symmetric/Hermitian;
5718:                             the relevant flags must be set independently.


5721:    Options For Use with MatSetValues():
5722:    Insert a logically dense subblock, which can be
5723: .    MAT_ROW_ORIENTED - row-oriented (default)

5725:    Note these options reflect the data you pass in with MatSetValues(); it has
5726:    nothing to do with how the data is stored internally in the matrix
5727:    data structure.

5729:    When (re)assembling a matrix, we can restrict the input for
5730:    efficiency/debugging purposes.  These options include:
5731: +    MAT_NEW_NONZERO_LOCATIONS - additional insertions will be allowed if they generate a new nonzero (slow)
5732: .    MAT_FORCE_DIAGONAL_ENTRIES - forces diagonal entries to be allocated
5733: .    MAT_IGNORE_OFF_PROC_ENTRIES - drops off-processor entries
5734: .    MAT_NEW_NONZERO_LOCATION_ERR - generates an error for new matrix entry
5735: .    MAT_USE_HASH_TABLE - uses a hash table to speed up matrix assembly
5736: .    MAT_NO_OFF_PROC_ENTRIES - you know each process will only set values for its own rows, will generate an error if
5737:         any process sets values for another process. This avoids all reductions in the MatAssembly routines and thus improves
5738:         performance for very large process counts.
5739: -    MAT_SUBSET_OFF_PROC_ENTRIES - you know that the first assembly after setting this flag will set a superset
5740:         of the off-process entries required for all subsequent assemblies. This avoids a rendezvous step in the MatAssembly
5741:         functions, instead sending only neighbor messages.

5743:    Notes:
5744:    Except for MAT_UNUSED_NONZERO_LOCATION_ERR and  MAT_ROW_ORIENTED all processes that share the matrix must pass the same value in flg!

5746:    Some options are relevant only for particular matrix types and
5747:    are thus ignored by others.  Other options are not supported by
5748:    certain matrix types and will generate an error message if set.

5750:    If using a Fortran 77 module to compute a matrix, one may need to
5751:    use the column-oriented option (or convert to the row-oriented
5752:    format).

5754:    MAT_NEW_NONZERO_LOCATIONS set to PETSC_FALSE indicates that any add or insertion
5755:    that would generate a new entry in the nonzero structure is instead
5756:    ignored.  Thus, if memory has not alredy been allocated for this particular
5757:    data, then the insertion is ignored. For dense matrices, in which
5758:    the entire array is allocated, no entries are ever ignored.
5759:    Set after the first MatAssemblyEnd(). If this option is set then the MatAssemblyBegin/End() processes has one less global reduction

5761:    MAT_NEW_NONZERO_LOCATION_ERR set to PETSC_TRUE indicates that any add or insertion
5762:    that would generate a new entry in the nonzero structure instead produces
5763:    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

5765:    MAT_NEW_NONZERO_ALLOCATION_ERR set to PETSC_TRUE indicates that any add or insertion
5766:    that would generate a new entry that has not been preallocated will
5767:    instead produce an error. (Currently supported for AIJ and BAIJ formats
5768:    only.) This is a useful flag when debugging matrix memory preallocation.
5769:    If this option is set then the MatAssemblyBegin/End() processes has one less global reduction

5771:    MAT_IGNORE_OFF_PROC_ENTRIES set to PETSC_TRUE indicates entries destined for
5772:    other processors should be dropped, rather than stashed.
5773:    This is useful if you know that the "owning" processor is also
5774:    always generating the correct matrix entries, so that PETSc need
5775:    not transfer duplicate entries generated on another processor.

5777:    MAT_USE_HASH_TABLE indicates that a hash table be used to improve the
5778:    searches during matrix assembly. When this flag is set, the hash table
5779:    is created during the first Matrix Assembly. This hash table is
5780:    used the next time through, during MatSetVaules()/MatSetVaulesBlocked()
5781:    to improve the searching of indices. MAT_NEW_NONZERO_LOCATIONS flag
5782:    should be used with MAT_USE_HASH_TABLE flag. This option is currently
5783:    supported by MATMPIBAIJ format only.

5785:    MAT_KEEP_NONZERO_PATTERN indicates when MatZeroRows() is called the zeroed entries
5786:    are kept in the nonzero structure

5788:    MAT_IGNORE_ZERO_ENTRIES - for AIJ/IS matrices this will stop zero values from creating
5789:    a zero location in the matrix

5791:    MAT_USE_INODES - indicates using inode version of the code - works with AIJ matrix types

5793:    MAT_NO_OFF_PROC_ZERO_ROWS - you know each process will only zero its own rows. This avoids all reductions in the
5794:         zero row routines and thus improves performance for very large process counts.

5796:    MAT_IGNORE_LOWER_TRIANGULAR - For SBAIJ matrices will ignore any insertions you make in the lower triangular
5797:         part of the matrix (since they should match the upper triangular part).

5799:    MAT_SORTED_FULL - each process provides exactly its local rows; all column indices for a given row are passed in a
5800:                      single call to MatSetValues(), preallocation is perfect, row oriented, INSERT_VALUES is used. Common
5801:                      with finite difference schemes with non-periodic boundary conditions.

5803:    Level: intermediate

5805: .seealso:  MatOption, Mat

5807: @*/
5808: PetscErrorCode MatSetOption(Mat mat,MatOption op,PetscBool flg)
5809: {

5814:   if (op > 0) {
5817:   }

5819:   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);

5821:   switch (op) {
5822:   case MAT_FORCE_DIAGONAL_ENTRIES:
5823:     mat->force_diagonals = flg;
5824:     return(0);
5825:   case MAT_NO_OFF_PROC_ENTRIES:
5826:     mat->nooffprocentries = flg;
5827:     return(0);
5828:   case MAT_SUBSET_OFF_PROC_ENTRIES:
5829:     mat->assembly_subset = flg;
5830:     if (!mat->assembly_subset) { /* See the same logic in VecAssembly wrt VEC_SUBSET_OFF_PROC_ENTRIES */
5831: #if !defined(PETSC_HAVE_MPIUNI)
5832:       MatStashScatterDestroy_BTS(&mat->stash);
5833: #endif
5834:       mat->stash.first_assembly_done = PETSC_FALSE;
5835:     }
5836:     return(0);
5837:   case MAT_NO_OFF_PROC_ZERO_ROWS:
5838:     mat->nooffproczerorows = flg;
5839:     return(0);
5840:   case MAT_SPD:
5841:     mat->spd_set = PETSC_TRUE;
5842:     mat->spd     = flg;
5843:     if (flg) {
5844:       mat->symmetric                  = PETSC_TRUE;
5845:       mat->structurally_symmetric     = PETSC_TRUE;
5846:       mat->symmetric_set              = PETSC_TRUE;
5847:       mat->structurally_symmetric_set = PETSC_TRUE;
5848:     }
5849:     break;
5850:   case MAT_SYMMETRIC:
5851:     mat->symmetric = flg;
5852:     if (flg) mat->structurally_symmetric = PETSC_TRUE;
5853:     mat->symmetric_set              = PETSC_TRUE;
5854:     mat->structurally_symmetric_set = flg;
5855: #if !defined(PETSC_USE_COMPLEX)
5856:     mat->hermitian     = flg;
5857:     mat->hermitian_set = PETSC_TRUE;
5858: #endif
5859:     break;
5860:   case MAT_HERMITIAN:
5861:     mat->hermitian = flg;
5862:     if (flg) mat->structurally_symmetric = PETSC_TRUE;
5863:     mat->hermitian_set              = PETSC_TRUE;
5864:     mat->structurally_symmetric_set = flg;
5865: #if !defined(PETSC_USE_COMPLEX)
5866:     mat->symmetric     = flg;
5867:     mat->symmetric_set = PETSC_TRUE;
5868: #endif
5869:     break;
5870:   case MAT_STRUCTURALLY_SYMMETRIC:
5871:     mat->structurally_symmetric     = flg;
5872:     mat->structurally_symmetric_set = PETSC_TRUE;
5873:     break;
5874:   case MAT_SYMMETRY_ETERNAL:
5875:     mat->symmetric_eternal = flg;
5876:     break;
5877:   case MAT_STRUCTURE_ONLY:
5878:     mat->structure_only = flg;
5879:     break;
5880:   case MAT_SORTED_FULL:
5881:     mat->sortedfull = flg;
5882:     break;
5883:   default:
5884:     break;
5885:   }
5886:   if (mat->ops->setoption) {
5887:     (*mat->ops->setoption)(mat,op,flg);
5888:   }
5889:   return(0);
5890: }

5892: /*@
5893:    MatGetOption - Gets a parameter option that has been set for a matrix.

5895:    Logically Collective on Mat for certain operations, such as MAT_SPD, not collective for MAT_ROW_ORIENTED, see MatOption

5897:    Input Parameters:
5898: +  mat - the matrix
5899: -  option - the option, this only responds to certain options, check the code for which ones

5901:    Output Parameter:
5902: .  flg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE)

5904:     Notes:
5905:     Can only be called after MatSetSizes() and MatSetType() have been set.

5907:    Level: intermediate

5909: .seealso:  MatOption, MatSetOption()

5911: @*/
5912: PetscErrorCode MatGetOption(Mat mat,MatOption op,PetscBool *flg)
5913: {

5918:   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);
5919:   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()");

5921:   switch (op) {
5922:   case MAT_NO_OFF_PROC_ENTRIES:
5923:     *flg = mat->nooffprocentries;
5924:     break;
5925:   case MAT_NO_OFF_PROC_ZERO_ROWS:
5926:     *flg = mat->nooffproczerorows;
5927:     break;
5928:   case MAT_SYMMETRIC:
5929:     *flg = mat->symmetric;
5930:     break;
5931:   case MAT_HERMITIAN:
5932:     *flg = mat->hermitian;
5933:     break;
5934:   case MAT_STRUCTURALLY_SYMMETRIC:
5935:     *flg = mat->structurally_symmetric;
5936:     break;
5937:   case MAT_SYMMETRY_ETERNAL:
5938:     *flg = mat->symmetric_eternal;
5939:     break;
5940:   case MAT_SPD:
5941:     *flg = mat->spd;
5942:     break;
5943:   default:
5944:     break;
5945:   }
5946:   return(0);
5947: }

5949: /*@
5950:    MatZeroEntries - Zeros all entries of a matrix.  For sparse matrices
5951:    this routine retains the old nonzero structure.

5953:    Logically Collective on Mat

5955:    Input Parameters:
5956: .  mat - the matrix

5958:    Level: intermediate

5960:    Notes:
5961:     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.
5962:    See the Performance chapter of the users manual for information on preallocating matrices.

5964: .seealso: MatZeroRows()
5965: @*/
5966: PetscErrorCode MatZeroEntries(Mat mat)
5967: {

5973:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
5974:   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");
5975:   if (!mat->ops->zeroentries) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
5976:   MatCheckPreallocated(mat,1);

5978:   PetscLogEventBegin(MAT_ZeroEntries,mat,0,0,0);
5979:   (*mat->ops->zeroentries)(mat);
5980:   PetscLogEventEnd(MAT_ZeroEntries,mat,0,0,0);
5981:   PetscObjectStateIncrease((PetscObject)mat);
5982:   return(0);
5983: }

5985: /*@
5986:    MatZeroRowsColumns - Zeros all entries (except possibly the main diagonal)
5987:    of a set of rows and columns of a matrix.

5989:    Collective on Mat

5991:    Input Parameters:
5992: +  mat - the matrix
5993: .  numRows - the number of rows to remove
5994: .  rows - the global row indices
5995: .  diag - value put in all diagonals of eliminated rows (0.0 will even eliminate diagonal entry)
5996: .  x - optional vector of solutions for zeroed rows (other entries in vector are not used)
5997: -  b - optional vector of right hand side, that will be adjusted by provided solution

5999:    Notes:
6000:    This does not change the nonzero structure of the matrix, it merely zeros those entries in the matrix.

6002:    The user can set a value in the diagonal entry (or for the AIJ and
6003:    row formats can optionally remove the main diagonal entry from the
6004:    nonzero structure as well, by passing 0.0 as the final argument).

6006:    For the parallel case, all processes that share the matrix (i.e.,
6007:    those in the communicator used for matrix creation) MUST call this
6008:    routine, regardless of whether any rows being zeroed are owned by
6009:    them.

6011:    Each processor can indicate any rows in the entire matrix to be zeroed (i.e. each process does NOT have to
6012:    list only rows local to itself).

6014:    The option MAT_NO_OFF_PROC_ZERO_ROWS does not apply to this routine.

6016:    Level: intermediate

6018: .seealso: MatZeroRowsIS(), MatZeroRows(), MatZeroRowsLocalIS(), MatZeroRowsStencil(), MatZeroEntries(), MatZeroRowsLocal(), MatSetOption(),
6019:           MatZeroRowsColumnsLocal(), MatZeroRowsColumnsLocalIS(), MatZeroRowsColumnsIS(), MatZeroRowsColumnsStencil()
6020: @*/
6021: PetscErrorCode MatZeroRowsColumns(Mat mat,PetscInt numRows,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
6022: {

6029:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
6030:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
6031:   if (!mat->ops->zerorowscolumns) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
6032:   MatCheckPreallocated(mat,1);

6034:   (*mat->ops->zerorowscolumns)(mat,numRows,rows,diag,x,b);
6035:   MatViewFromOptions(mat,NULL,"-mat_view");
6036:   PetscObjectStateIncrease((PetscObject)mat);
6037:   return(0);
6038: }

6040: /*@
6041:    MatZeroRowsColumnsIS - Zeros all entries (except possibly the main diagonal)
6042:    of a set of rows and columns of a matrix.

6044:    Collective on Mat

6046:    Input Parameters:
6047: +  mat - the matrix
6048: .  is - the rows to zero
6049: .  diag - value put in all diagonals of eliminated rows (0.0 will even eliminate diagonal entry)
6050: .  x - optional vector of solutions for zeroed rows (other entries in vector are not used)
6051: -  b - optional vector of right hand side, that will be adjusted by provided solution

6053:    Notes:
6054:    This does not change the nonzero structure of the matrix, it merely zeros those entries in the matrix.

6056:    The user can set a value in the diagonal entry (or for the AIJ and
6057:    row formats can optionally remove the main diagonal entry from the
6058:    nonzero structure as well, by passing 0.0 as the final argument).

6060:    For the parallel case, all processes that share the matrix (i.e.,
6061:    those in the communicator used for matrix creation) MUST call this
6062:    routine, regardless of whether any rows being zeroed are owned by
6063:    them.

6065:    Each processor can indicate any rows in the entire matrix to be zeroed (i.e. each process does NOT have to
6066:    list only rows local to itself).

6068:    The option MAT_NO_OFF_PROC_ZERO_ROWS does not apply to this routine.

6070:    Level: intermediate

6072: .seealso: MatZeroRowsIS(), MatZeroRowsColumns(), MatZeroRowsLocalIS(), MatZeroRowsStencil(), MatZeroEntries(), MatZeroRowsLocal(), MatSetOption(),
6073:           MatZeroRowsColumnsLocal(), MatZeroRowsColumnsLocalIS(), MatZeroRows(), MatZeroRowsColumnsStencil()
6074: @*/
6075: PetscErrorCode MatZeroRowsColumnsIS(Mat mat,IS is,PetscScalar diag,Vec x,Vec b)
6076: {
6078:   PetscInt       numRows;
6079:   const PetscInt *rows;

6086:   ISGetLocalSize(is,&numRows);
6087:   ISGetIndices(is,&rows);
6088:   MatZeroRowsColumns(mat,numRows,rows,diag,x,b);
6089:   ISRestoreIndices(is,&rows);
6090:   return(0);
6091: }

6093: /*@
6094:    MatZeroRows - Zeros all entries (except possibly the main diagonal)
6095:    of a set of rows of a matrix.

6097:    Collective on Mat

6099:    Input Parameters:
6100: +  mat - the matrix
6101: .  numRows - the number of rows to remove
6102: .  rows - the global row indices
6103: .  diag - value put in all diagonals of eliminated rows (0.0 will even eliminate diagonal entry)
6104: .  x - optional vector of solutions for zeroed rows (other entries in vector are not used)
6105: -  b - optional vector of right hand side, that will be adjusted by provided solution

6107:    Notes:
6108:    For the AIJ and BAIJ matrix formats this removes the old nonzero structure,
6109:    but does not release memory.  For the dense and block diagonal
6110:    formats this does not alter the nonzero structure.

6112:    If the option MatSetOption(mat,MAT_KEEP_NONZERO_PATTERN,PETSC_TRUE) the nonzero structure
6113:    of the matrix is not changed (even for AIJ and BAIJ matrices) the values are
6114:    merely zeroed.

6116:    The user can set a value in the diagonal entry (or for the AIJ and
6117:    row formats can optionally remove the main diagonal entry from the
6118:    nonzero structure as well, by passing 0.0 as the final argument).

6120:    For the parallel case, all processes that share the matrix (i.e.,
6121:    those in the communicator used for matrix creation) MUST call this
6122:    routine, regardless of whether any rows being zeroed are owned by
6123:    them.

6125:    Each processor can indicate any rows in the entire matrix to be zeroed (i.e. each process does NOT have to
6126:    list only rows local to itself).

6128:    You can call MatSetOption(mat,MAT_NO_OFF_PROC_ZERO_ROWS,PETSC_TRUE) if each process indicates only rows it
6129:    owns that are to be zeroed. This saves a global synchronization in the implementation.

6131:    Level: intermediate

6133: .seealso: MatZeroRowsIS(), MatZeroRowsColumns(), MatZeroRowsLocalIS(), MatZeroRowsStencil(), MatZeroEntries(), MatZeroRowsLocal(), MatSetOption(),
6134:           MatZeroRowsColumnsLocal(), MatZeroRowsColumnsLocalIS(), MatZeroRowsColumnsIS(), MatZeroRowsColumnsStencil()
6135: @*/
6136: PetscErrorCode MatZeroRows(Mat mat,PetscInt numRows,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
6137: {

6144:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
6145:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
6146:   if (!mat->ops->zerorows) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
6147:   MatCheckPreallocated(mat,1);

6149:   (*mat->ops->zerorows)(mat,numRows,rows,diag,x,b);
6150:   MatViewFromOptions(mat,NULL,"-mat_view");
6151:   PetscObjectStateIncrease((PetscObject)mat);
6152:   return(0);
6153: }

6155: /*@
6156:    MatZeroRowsIS - Zeros all entries (except possibly the main diagonal)
6157:    of a set of rows of a matrix.

6159:    Collective on Mat

6161:    Input Parameters:
6162: +  mat - the matrix
6163: .  is - index set of rows to remove (if NULL then no row is removed)
6164: .  diag - value put in all diagonals of eliminated rows
6165: .  x - optional vector of solutions for zeroed rows (other entries in vector are not used)
6166: -  b - optional vector of right hand side, that will be adjusted by provided solution

6168:    Notes:
6169:    For the AIJ and BAIJ matrix formats this removes the old nonzero structure,
6170:    but does not release memory.  For the dense and block diagonal
6171:    formats this does not alter the nonzero structure.

6173:    If the option MatSetOption(mat,MAT_KEEP_NONZERO_PATTERN,PETSC_TRUE) the nonzero structure
6174:    of the matrix is not changed (even for AIJ and BAIJ matrices) the values are
6175:    merely zeroed.

6177:    The user can set a value in the diagonal entry (or for the AIJ and
6178:    row formats can optionally remove the main diagonal entry from the
6179:    nonzero structure as well, by passing 0.0 as the final argument).

6181:    For the parallel case, all processes that share the matrix (i.e.,
6182:    those in the communicator used for matrix creation) MUST call this
6183:    routine, regardless of whether any rows being zeroed are owned by
6184:    them.

6186:    Each processor can indicate any rows in the entire matrix to be zeroed (i.e. each process does NOT have to
6187:    list only rows local to itself).

6189:    You can call MatSetOption(mat,MAT_NO_OFF_PROC_ZERO_ROWS,PETSC_TRUE) if each process indicates only rows it
6190:    owns that are to be zeroed. This saves a global synchronization in the implementation.

6192:    Level: intermediate

6194: .seealso: MatZeroRows(), MatZeroRowsColumns(), MatZeroRowsLocalIS(), MatZeroRowsStencil(), MatZeroEntries(), MatZeroRowsLocal(), MatSetOption(),
6195:           MatZeroRowsColumnsLocal(), MatZeroRowsColumnsLocalIS(), MatZeroRowsColumnsIS(), MatZeroRowsColumnsStencil()
6196: @*/
6197: PetscErrorCode MatZeroRowsIS(Mat mat,IS is,PetscScalar diag,Vec x,Vec b)
6198: {
6199:   PetscInt       numRows = 0;
6200:   const PetscInt *rows = NULL;

6206:   if (is) {
6208:     ISGetLocalSize(is,&numRows);
6209:     ISGetIndices(is,&rows);
6210:   }
6211:   MatZeroRows(mat,numRows,rows,diag,x,b);
6212:   if (is) {
6213:     ISRestoreIndices(is,&rows);
6214:   }
6215:   return(0);
6216: }

6218: /*@
6219:    MatZeroRowsStencil - Zeros all entries (except possibly the main diagonal)
6220:    of a set of rows of a matrix. These rows must be local to the process.

6222:    Collective on Mat

6224:    Input Parameters:
6225: +  mat - the matrix
6226: .  numRows - the number of rows to remove
6227: .  rows - the grid coordinates (and component number when dof > 1) for matrix rows
6228: .  diag - value put in all diagonals of eliminated rows (0.0 will even eliminate diagonal entry)
6229: .  x - optional vector of solutions for zeroed rows (other entries in vector are not used)
6230: -  b - optional vector of right hand side, that will be adjusted by provided solution

6232:    Notes:
6233:    For the AIJ and BAIJ matrix formats this removes the old nonzero structure,
6234:    but does not release memory.  For the dense and block diagonal
6235:    formats this does not alter the nonzero structure.

6237:    If the option MatSetOption(mat,MAT_KEEP_NONZERO_PATTERN,PETSC_TRUE) the nonzero structure
6238:    of the matrix is not changed (even for AIJ and BAIJ matrices) the values are
6239:    merely zeroed.

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:    For the parallel case, all processes that share the matrix (i.e.,
6246:    those in the communicator used for matrix creation) MUST call this
6247:    routine, regardless of whether any rows being zeroed are owned by
6248:    them.

6250:    Each processor can indicate any rows in the entire matrix to be zeroed (i.e. each process does NOT have to
6251:    list only rows local to itself).

6253:    The grid coordinates are across the entire grid, not just the local portion

6255:    In Fortran idxm and idxn should be declared as
6256: $     MatStencil idxm(4,m)
6257:    and the values inserted using
6258: $    idxm(MatStencil_i,1) = i
6259: $    idxm(MatStencil_j,1) = j
6260: $    idxm(MatStencil_k,1) = k
6261: $    idxm(MatStencil_c,1) = c
6262:    etc

6264:    For periodic boundary conditions use negative indices for values to the left (below 0; that are to be
6265:    obtained by wrapping values from right edge). For values to the right of the last entry using that index plus one
6266:    etc to obtain values that obtained by wrapping the values from the left edge. This does not work for anything but the
6267:    DM_BOUNDARY_PERIODIC boundary type.

6269:    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
6270:    a single value per point) you can skip filling those indices.

6272:    Level: intermediate

6274: .seealso: MatZeroRowsIS(), MatZeroRowsColumns(), MatZeroRowsLocalIS(), MatZeroRowsl(), MatZeroEntries(), MatZeroRowsLocal(), MatSetOption(),
6275:           MatZeroRowsColumnsLocal(), MatZeroRowsColumnsLocalIS(), MatZeroRowsColumnsIS(), MatZeroRowsColumnsStencil()
6276: @*/
6277: PetscErrorCode MatZeroRowsStencil(Mat mat,PetscInt numRows,const MatStencil rows[],PetscScalar diag,Vec x,Vec b)
6278: {
6279:   PetscInt       dim     = mat->stencil.dim;
6280:   PetscInt       sdim    = dim - (1 - (PetscInt) mat->stencil.noc);
6281:   PetscInt       *dims   = mat->stencil.dims+1;
6282:   PetscInt       *starts = mat->stencil.starts;
6283:   PetscInt       *dxm    = (PetscInt*) rows;
6284:   PetscInt       *jdxm, i, j, tmp, numNewRows = 0;


6292:   PetscMalloc1(numRows, &jdxm);
6293:   for (i = 0; i < numRows; ++i) {
6294:     /* Skip unused dimensions (they are ordered k, j, i, c) */
6295:     for (j = 0; j < 3-sdim; ++j) dxm++;
6296:     /* Local index in X dir */
6297:     tmp = *dxm++ - starts[0];
6298:     /* Loop over remaining dimensions */
6299:     for (j = 0; j < dim-1; ++j) {
6300:       /* If nonlocal, set index to be negative */
6301:       if ((*dxm++ - starts[j+1]) < 0 || tmp < 0) tmp = PETSC_MIN_INT;
6302:       /* Update local index */
6303:       else tmp = tmp*dims[j] + *(dxm-1) - starts[j+1];
6304:     }
6305:     /* Skip component slot if necessary */
6306:     if (mat->stencil.noc) dxm++;
6307:     /* Local row number */
6308:     if (tmp >= 0) {
6309:       jdxm[numNewRows++] = tmp;
6310:     }
6311:   }
6312:   MatZeroRowsLocal(mat,numNewRows,jdxm,diag,x,b);
6313:   PetscFree(jdxm);
6314:   return(0);
6315: }

6317: /*@
6318:    MatZeroRowsColumnsStencil - Zeros all row and column entries (except possibly the main diagonal)
6319:    of a set of rows and columns of a matrix.

6321:    Collective on Mat

6323:    Input Parameters:
6324: +  mat - the matrix
6325: .  numRows - the number of rows/columns to remove
6326: .  rows - the grid coordinates (and component number when dof > 1) for matrix rows
6327: .  diag - value put in all diagonals of eliminated rows (0.0 will even eliminate diagonal entry)
6328: .  x - optional vector of solutions for zeroed rows (other entries in vector are not used)
6329: -  b - optional vector of right hand side, that will be adjusted by provided solution

6331:    Notes:
6332:    For the AIJ and BAIJ matrix formats this removes the old nonzero structure,
6333:    but does not release memory.  For the dense and block diagonal
6334:    formats this does not alter the nonzero structure.

6336:    If the option MatSetOption(mat,MAT_KEEP_NONZERO_PATTERN,PETSC_TRUE) the nonzero structure
6337:    of the matrix is not changed (even for AIJ and BAIJ matrices) the values are
6338:    merely zeroed.

6340:    The user can set a value in the diagonal entry (or for the AIJ and
6341:    row formats can optionally remove the main diagonal entry from the
6342:    nonzero structure as well, by passing 0.0 as the final argument).

6344:    For the parallel case, all processes that share the matrix (i.e.,
6345:    those in the communicator used for matrix creation) MUST call this
6346:    routine, regardless of whether any rows being zeroed are owned by
6347:    them.

6349:    Each processor can indicate any rows in the entire matrix to be zeroed (i.e. each process does NOT have to
6350:    list only rows local to itself, but the row/column numbers are given in local numbering).

6352:    The grid coordinates are across the entire grid, not just the local portion

6354:    In Fortran idxm and idxn should be declared as
6355: $     MatStencil idxm(4,m)
6356:    and the values inserted using
6357: $    idxm(MatStencil_i,1) = i
6358: $    idxm(MatStencil_j,1) = j
6359: $    idxm(MatStencil_k,1) = k
6360: $    idxm(MatStencil_c,1) = c
6361:    etc

6363:    For periodic boundary conditions use negative indices for values to the left (below 0; that are to be
6364:    obtained by wrapping values from right edge). For values to the right of the last entry using that index plus one
6365:    etc to obtain values that obtained by wrapping the values from the left edge. This does not work for anything but the
6366:    DM_BOUNDARY_PERIODIC boundary type.

6368:    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
6369:    a single value per point) you can skip filling those indices.

6371:    Level: intermediate

6373: .seealso: MatZeroRowsIS(), MatZeroRowsColumns(), MatZeroRowsLocalIS(), MatZeroRowsStencil(), MatZeroEntries(), MatZeroRowsLocal(), MatSetOption(),
6374:           MatZeroRowsColumnsLocal(), MatZeroRowsColumnsLocalIS(), MatZeroRowsColumnsIS(), MatZeroRows()
6375: @*/
6376: PetscErrorCode MatZeroRowsColumnsStencil(Mat mat,PetscInt numRows,const MatStencil rows[],PetscScalar diag,Vec x,Vec b)
6377: {
6378:   PetscInt       dim     = mat->stencil.dim;
6379:   PetscInt       sdim    = dim - (1 - (PetscInt) mat->stencil.noc);
6380:   PetscInt       *dims   = mat->stencil.dims+1;
6381:   PetscInt       *starts = mat->stencil.starts;
6382:   PetscInt       *dxm    = (PetscInt*) rows;
6383:   PetscInt       *jdxm, i, j, tmp, numNewRows = 0;


6391:   PetscMalloc1(numRows, &jdxm);
6392:   for (i = 0; i < numRows; ++i) {
6393:     /* Skip unused dimensions (they are ordered k, j, i, c) */
6394:     for (j = 0; j < 3-sdim; ++j) dxm++;
6395:     /* Local index in X dir */
6396:     tmp = *dxm++ - starts[0];
6397:     /* Loop over remaining dimensions */
6398:     for (j = 0; j < dim-1; ++j) {
6399:       /* If nonlocal, set index to be negative */
6400:       if ((*dxm++ - starts[j+1]) < 0 || tmp < 0) tmp = PETSC_MIN_INT;
6401:       /* Update local index */
6402:       else tmp = tmp*dims[j] + *(dxm-1) - starts[j+1];
6403:     }
6404:     /* Skip component slot if necessary */
6405:     if (mat->stencil.noc) dxm++;
6406:     /* Local row number */
6407:     if (tmp >= 0) {
6408:       jdxm[numNewRows++] = tmp;
6409:     }
6410:   }
6411:   MatZeroRowsColumnsLocal(mat,numNewRows,jdxm,diag,x,b);
6412:   PetscFree(jdxm);
6413:   return(0);
6414: }

6416: /*@C
6417:    MatZeroRowsLocal - Zeros all entries (except possibly the main diagonal)
6418:    of a set of rows of a matrix; using local numbering of rows.

6420:    Collective on Mat

6422:    Input Parameters:
6423: +  mat - the matrix
6424: .  numRows - the number of rows to remove
6425: .  rows - the local row indices
6426: .  diag - value put in all diagonals of eliminated rows
6427: .  x - optional vector of solutions for zeroed rows (other entries in vector are not used)
6428: -  b - optional vector of right hand side, that will be adjusted by provided solution

6430:    Notes:
6431:    Before calling MatZeroRowsLocal(), the user must first set the
6432:    local-to-global mapping by calling MatSetLocalToGlobalMapping().

6434:    For the AIJ matrix formats this removes the old nonzero structure,
6435:    but does not release memory.  For the dense and block diagonal
6436:    formats this does not alter the nonzero structure.

6438:    If the option MatSetOption(mat,MAT_KEEP_NONZERO_PATTERN,PETSC_TRUE) the nonzero structure
6439:    of the matrix is not changed (even for AIJ and BAIJ matrices) the values are
6440:    merely zeroed.

6442:    The user can set a value in the diagonal entry (or for the AIJ and
6443:    row formats can optionally remove the main diagonal entry from the
6444:    nonzero structure as well, by passing 0.0 as the final argument).

6446:    You can call MatSetOption(mat,MAT_NO_OFF_PROC_ZERO_ROWS,PETSC_TRUE) if each process indicates only rows it
6447:    owns that are to be zeroed. This saves a global synchronization in the implementation.

6449:    Level: intermediate

6451: .seealso: MatZeroRowsIS(), MatZeroRowsColumns(), MatZeroRowsLocalIS(), MatZeroRowsStencil(), MatZeroEntries(), MatZeroRows(), MatSetOption(),
6452:           MatZeroRowsColumnsLocal(), MatZeroRowsColumnsLocalIS(), MatZeroRowsColumnsIS(), MatZeroRowsColumnsStencil()
6453: @*/
6454: PetscErrorCode MatZeroRowsLocal(Mat mat,PetscInt numRows,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
6455: {

6462:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
6463:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
6464:   MatCheckPreallocated(mat,1);

6466:   if (mat->ops->zerorowslocal) {
6467:     (*mat->ops->zerorowslocal)(mat,numRows,rows,diag,x,b);
6468:   } else {
6469:     IS             is, newis;
6470:     const PetscInt *newRows;

6472:     if (!mat->rmap->mapping) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Need to provide local to global mapping to matrix first");
6473:     ISCreateGeneral(PETSC_COMM_SELF,numRows,rows,PETSC_COPY_VALUES,&is);
6474:     ISLocalToGlobalMappingApplyIS(mat->rmap->mapping,is,&newis);
6475:     ISGetIndices(newis,&newRows);
6476:     (*mat->ops->zerorows)(mat,numRows,newRows,diag,x,b);
6477:     ISRestoreIndices(newis,&newRows);
6478:     ISDestroy(&newis);
6479:     ISDestroy(&is);
6480:   }
6481:   PetscObjectStateIncrease((PetscObject)mat);
6482:   return(0);
6483: }

6485: /*@
6486:    MatZeroRowsLocalIS - Zeros all entries (except possibly the main diagonal)
6487:    of a set of rows of a matrix; using local numbering of rows.

6489:    Collective on Mat

6491:    Input Parameters:
6492: +  mat - the matrix
6493: .  is - index set of rows to remove
6494: .  diag - value put in all diagonals of eliminated rows
6495: .  x - optional vector of solutions for zeroed rows (other entries in vector are not used)
6496: -  b - optional vector of right hand side, that will be adjusted by provided solution

6498:    Notes:
6499:    Before calling MatZeroRowsLocalIS(), the user must first set the
6500:    local-to-global mapping by calling MatSetLocalToGlobalMapping().

6502:    For the AIJ matrix formats this removes the old nonzero structure,
6503:    but does not release memory.  For the dense and block diagonal
6504:    formats this does not alter the nonzero structure.

6506:    If the option MatSetOption(mat,MAT_KEEP_NONZERO_PATTERN,PETSC_TRUE) the nonzero structure
6507:    of the matrix is not changed (even for AIJ and BAIJ matrices) the values are
6508:    merely zeroed.

6510:    The user can set a value in the diagonal entry (or for the AIJ and
6511:    row formats can optionally remove the main diagonal entry from the
6512:    nonzero structure as well, by passing 0.0 as the final argument).

6514:    You can call MatSetOption(mat,MAT_NO_OFF_PROC_ZERO_ROWS,PETSC_TRUE) if each process indicates only rows it
6515:    owns that are to be zeroed. This saves a global synchronization in the implementation.

6517:    Level: intermediate

6519: .seealso: MatZeroRowsIS(), MatZeroRowsColumns(), MatZeroRows(), MatZeroRowsStencil(), MatZeroEntries(), MatZeroRowsLocal(), MatSetOption(),
6520:           MatZeroRowsColumnsLocal(), MatZeroRowsColumnsLocalIS(), MatZeroRowsColumnsIS(), MatZeroRowsColumnsStencil()
6521: @*/
6522: PetscErrorCode MatZeroRowsLocalIS(Mat mat,IS is,PetscScalar diag,Vec x,Vec b)
6523: {
6525:   PetscInt       numRows;
6526:   const PetscInt *rows;

6532:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
6533:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
6534:   MatCheckPreallocated(mat,1);

6536:   ISGetLocalSize(is,&numRows);
6537:   ISGetIndices(is,&rows);
6538:   MatZeroRowsLocal(mat,numRows,rows,diag,x,b);
6539:   ISRestoreIndices(is,&rows);
6540:   return(0);
6541: }

6543: /*@
6544:    MatZeroRowsColumnsLocal - Zeros all entries (except possibly the main diagonal)
6545:    of a set of rows and columns of a matrix; using local numbering of rows.

6547:    Collective on Mat

6549:    Input Parameters:
6550: +  mat - the matrix
6551: .  numRows - the number of rows to remove
6552: .  rows - the global row indices
6553: .  diag - value put in all diagonals of eliminated rows
6554: .  x - optional vector of solutions for zeroed rows (other entries in vector are not used)
6555: -  b - optional vector of right hand side, that will be adjusted by provided solution

6557:    Notes:
6558:    Before calling MatZeroRowsColumnsLocal(), the user must first set the
6559:    local-to-global mapping by calling MatSetLocalToGlobalMapping().

6561:    The user can set a value in the diagonal entry (or for the AIJ and
6562:    row formats can optionally remove the main diagonal entry from the
6563:    nonzero structure as well, by passing 0.0 as the final argument).

6565:    Level: intermediate

6567: .seealso: MatZeroRowsIS(), MatZeroRowsColumns(), MatZeroRowsLocalIS(), MatZeroRowsStencil(), MatZeroEntries(), MatZeroRowsLocal(), MatSetOption(),
6568:           MatZeroRows(), MatZeroRowsColumnsLocalIS(), MatZeroRowsColumnsIS(), MatZeroRowsColumnsStencil()
6569: @*/
6570: PetscErrorCode MatZeroRowsColumnsLocal(Mat mat,PetscInt numRows,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
6571: {
6573:   IS             is, newis;
6574:   const PetscInt *newRows;

6580:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
6581:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
6582:   MatCheckPreallocated(mat,1);

6584:   if (!mat->cmap->mapping) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Need to provide local to global mapping to matrix first");
6585:   ISCreateGeneral(PETSC_COMM_SELF,numRows,rows,PETSC_COPY_VALUES,&is);
6586:   ISLocalToGlobalMappingApplyIS(mat->cmap->mapping,is,&newis);
6587:   ISGetIndices(newis,&newRows);
6588:   (*mat->ops->zerorowscolumns)(mat,numRows,newRows,diag,x,b);
6589:   ISRestoreIndices(newis,&newRows);
6590:   ISDestroy(&newis);
6591:   ISDestroy(&is);
6592:   PetscObjectStateIncrease((PetscObject)mat);
6593:   return(0);
6594: }

6596: /*@
6597:    MatZeroRowsColumnsLocalIS - Zeros all entries (except possibly the main diagonal)
6598:    of a set of rows and columns of a matrix; using local numbering of rows.

6600:    Collective on Mat

6602:    Input Parameters:
6603: +  mat - the matrix
6604: .  is - index set of rows to remove
6605: .  diag - value put in all diagonals of eliminated rows
6606: .  x - optional vector of solutions for zeroed rows (other entries in vector are not used)
6607: -  b - optional vector of right hand side, that will be adjusted by provided solution

6609:    Notes:
6610:    Before calling MatZeroRowsColumnsLocalIS(), the user must first set the
6611:    local-to-global mapping by calling MatSetLocalToGlobalMapping().

6613:    The user can set a value in the diagonal entry (or for the AIJ and
6614:    row formats can optionally remove the main diagonal entry from the
6615:    nonzero structure as well, by passing 0.0 as the final argument).

6617:    Level: intermediate

6619: .seealso: MatZeroRowsIS(), MatZeroRowsColumns(), MatZeroRowsLocalIS(), MatZeroRowsStencil(), MatZeroEntries(), MatZeroRowsLocal(), MatSetOption(),
6620:           MatZeroRowsColumnsLocal(), MatZeroRows(), MatZeroRowsColumnsIS(), MatZeroRowsColumnsStencil()
6621: @*/
6622: PetscErrorCode MatZeroRowsColumnsLocalIS(Mat mat,IS is,PetscScalar diag,Vec x,Vec b)
6623: {
6625:   PetscInt       numRows;
6626:   const PetscInt *rows;

6632:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
6633:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
6634:   MatCheckPreallocated(mat,1);

6636:   ISGetLocalSize(is,&numRows);
6637:   ISGetIndices(is,&rows);
6638:   MatZeroRowsColumnsLocal(mat,numRows,rows,diag,x,b);
6639:   ISRestoreIndices(is,&rows);
6640:   return(0);
6641: }

6643: /*@C
6644:    MatGetSize - Returns the numbers of rows and columns in a matrix.

6646:    Not Collective

6648:    Input Parameter:
6649: .  mat - the matrix

6651:    Output Parameters:
6652: +  m - the number of global rows
6653: -  n - the number of global columns

6655:    Note: both output parameters can be NULL on input.

6657:    Level: beginner

6659: .seealso: MatGetLocalSize()
6660: @*/
6661: PetscErrorCode MatGetSize(Mat mat,PetscInt *m,PetscInt *n)
6662: {
6665:   if (m) *m = mat->rmap->N;
6666:   if (n) *n = mat->cmap->N;
6667:   return(0);
6668: }

6670: /*@C
6671:    MatGetLocalSize - Returns the number of local rows and local columns
6672:    of a matrix, that is the local size of the left and right vectors as returned by MatCreateVecs().

6674:    Not Collective

6676:    Input Parameters:
6677: .  mat - the matrix

6679:    Output Parameters:
6680: +  m - the number of local rows
6681: -  n - the number of local columns

6683:    Note: both output parameters can be NULL on input.

6685:    Level: beginner

6687: .seealso: MatGetSize()
6688: @*/
6689: PetscErrorCode MatGetLocalSize(Mat mat,PetscInt *m,PetscInt *n)
6690: {
6695:   if (m) *m = mat->rmap->n;
6696:   if (n) *n = mat->cmap->n;
6697:   return(0);
6698: }

6700: /*@C
6701:    MatGetOwnershipRangeColumn - Returns the range of matrix columns associated with rows of a vector one multiplies by that owned by
6702:    this processor. (The columns of the "diagonal block")

6704:    Not Collective, unless matrix has not been allocated, then collective on Mat

6706:    Input Parameters:
6707: .  mat - the matrix

6709:    Output Parameters:
6710: +  m - the global index of the first local column
6711: -  n - one more than the global index of the last local column

6713:    Notes:
6714:     both output parameters can be NULL on input.

6716:    Level: developer

6718: .seealso:  MatGetOwnershipRange(), MatGetOwnershipRanges(), MatGetOwnershipRangesColumn()

6720: @*/
6721: PetscErrorCode MatGetOwnershipRangeColumn(Mat mat,PetscInt *m,PetscInt *n)
6722: {
6728:   MatCheckPreallocated(mat,1);
6729:   if (m) *m = mat->cmap->rstart;
6730:   if (n) *n = mat->cmap->rend;
6731:   return(0);
6732: }

6734: /*@C
6735:    MatGetOwnershipRange - Returns the range of matrix rows owned by
6736:    this processor, assuming that the matrix is laid out with the first
6737:    n1 rows on the first processor, the next n2 rows on the second, etc.
6738:    For certain parallel layouts this range may not be well defined.

6740:    Not Collective

6742:    Input Parameters:
6743: .  mat - the matrix

6745:    Output Parameters:
6746: +  m - the global index of the first local row
6747: -  n - one more than the global index of the last local row

6749:    Note: Both output parameters can be NULL on input.
6750: $  This function requires that the matrix be preallocated. If you have not preallocated, consider using
6751: $    PetscSplitOwnership(MPI_Comm comm, PetscInt *n, PetscInt *N)
6752: $  and then MPI_Scan() to calculate prefix sums of the local sizes.

6754:    Level: beginner

6756: .seealso:   MatGetOwnershipRanges(), MatGetOwnershipRangeColumn(), MatGetOwnershipRangesColumn(), PetscSplitOwnership(), PetscSplitOwnershipBlock()

6758: @*/
6759: PetscErrorCode MatGetOwnershipRange(Mat mat,PetscInt *m,PetscInt *n)
6760: {
6766:   MatCheckPreallocated(mat,1);
6767:   if (m) *m = mat->rmap->rstart;
6768:   if (n) *n = mat->rmap->rend;
6769:   return(0);
6770: }

6772: /*@C
6773:    MatGetOwnershipRanges - Returns the range of matrix rows owned by
6774:    each process

6776:    Not Collective, unless matrix has not been allocated, then collective on Mat

6778:    Input Parameters:
6779: .  mat - the matrix

6781:    Output Parameters:
6782: .  ranges - start of each processors portion plus one more than the total length at the end

6784:    Level: beginner

6786: .seealso:   MatGetOwnershipRange(), MatGetOwnershipRangeColumn(), MatGetOwnershipRangesColumn()

6788: @*/
6789: PetscErrorCode MatGetOwnershipRanges(Mat mat,const PetscInt **ranges)
6790: {

6796:   MatCheckPreallocated(mat,1);
6797:   PetscLayoutGetRanges(mat->rmap,ranges);
6798:   return(0);
6799: }

6801: /*@C
6802:    MatGetOwnershipRangesColumn - Returns the range of matrix columns associated with rows of a vector one multiplies by that owned by
6803:    this processor. (The columns of the "diagonal blocks" for each process)

6805:    Not Collective, unless matrix has not been allocated, then collective on Mat

6807:    Input Parameters:
6808: .  mat - the matrix

6810:    Output Parameters:
6811: .  ranges - start of each processors portion plus one more then the total length at the end

6813:    Level: beginner

6815: .seealso:   MatGetOwnershipRange(), MatGetOwnershipRangeColumn(), MatGetOwnershipRanges()

6817: @*/
6818: PetscErrorCode MatGetOwnershipRangesColumn(Mat mat,const PetscInt **ranges)
6819: {

6825:   MatCheckPreallocated(mat,1);
6826:   PetscLayoutGetRanges(mat->cmap,ranges);
6827:   return(0);
6828: }

6830: /*@C
6831:    MatGetOwnershipIS - Get row and column ownership as index sets

6833:    Not Collective

6835:    Input Arguments:
6836: .  A - matrix of type Elemental or ScaLAPACK

6838:    Output Arguments:
6839: +  rows - rows in which this process owns elements
6840: -  cols - columns in which this process owns elements

6842:    Level: intermediate

6844: .seealso: MatGetOwnershipRange(), MatGetOwnershipRangeColumn(), MatSetValues(), MATELEMENTAL
6845: @*/
6846: PetscErrorCode MatGetOwnershipIS(Mat A,IS *rows,IS *cols)
6847: {
6848:   PetscErrorCode ierr,(*f)(Mat,IS*,IS*);

6851:   MatCheckPreallocated(A,1);
6852:   PetscObjectQueryFunction((PetscObject)A,"MatGetOwnershipIS_C",&f);
6853:   if (f) {
6854:     (*f)(A,rows,cols);
6855:   } else {   /* Create a standard row-based partition, each process is responsible for ALL columns in their row block */
6856:     if (rows) {ISCreateStride(PETSC_COMM_SELF,A->rmap->n,A->rmap->rstart,1,rows);}
6857:     if (cols) {ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,cols);}
6858:   }
6859:   return(0);
6860: }

6862: /*@C
6863:    MatILUFactorSymbolic - Performs symbolic ILU factorization of a matrix.
6864:    Uses levels of fill only, not drop tolerance. Use MatLUFactorNumeric()
6865:    to complete the factorization.

6867:    Collective on Mat

6869:    Input Parameters:
6870: +  mat - the matrix
6871: .  row - row permutation
6872: .  column - column permutation
6873: -  info - structure containing
6874: $      levels - number of levels of fill.
6875: $      expected fill - as ratio of original fill.
6876: $      1 or 0 - indicating force fill on diagonal (improves robustness for matrices
6877:                 missing diagonal entries)

6879:    Output Parameters:
6880: .  fact - new matrix that has been symbolically factored

6882:    Notes:
6883:     See Users-Manual: ch_mat for additional information about choosing the fill factor for better efficiency.

6885:    Most users should employ the simplified KSP interface for linear solvers
6886:    instead of working directly with matrix algebra routines such as this.
6887:    See, e.g., KSPCreate().

6889:    Level: developer

6891: .seealso: MatLUFactorSymbolic(), MatLUFactorNumeric(), MatCholeskyFactor()
6892:           MatGetOrdering(), MatFactorInfo

6894:     Note: this uses the definition of level of fill as in Y. Saad, 2003

6896:     Developer Note: fortran interface is not autogenerated as the f90
6897:     interface defintion cannot be generated correctly [due to MatFactorInfo]

6899:    References:
6900:      Y. Saad, Iterative methods for sparse linear systems Philadelphia: Society for Industrial and Applied Mathematics, 2003
6901: @*/
6902: PetscErrorCode MatILUFactorSymbolic(Mat fact,Mat mat,IS row,IS col,const MatFactorInfo *info)
6903: {

6913:   if (info->levels < 0) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_OUTOFRANGE,"Levels of fill negative %D",(PetscInt)info->levels);
6914:   if (info->fill < 1.0) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_OUTOFRANGE,"Expected fill less than 1.0 %g",(double)info->fill);
6915:   if (!fact->ops->ilufactorsymbolic) {
6916:     MatSolverType stype;
6917:     MatFactorGetSolverType(fact,&stype);
6918:     SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Matrix type %s symbolic ILU using solver type %s",((PetscObject)mat)->type_name,stype);
6919:   }
6920:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
6921:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
6922:   MatCheckPreallocated(mat,2);

6924:   PetscLogEventBegin(MAT_ILUFactorSymbolic,mat,row,col,0);
6925:   (fact->ops->ilufactorsymbolic)(fact,mat,row,col,info);
6926:   PetscLogEventEnd(MAT_ILUFactorSymbolic,mat,row,col,0);
6927:   return(0);
6928: }

6930: /*@C
6931:    MatICCFactorSymbolic - Performs symbolic incomplete
6932:    Cholesky factorization for a symmetric matrix.  Use
6933:    MatCholeskyFactorNumeric() to complete the factorization.

6935:    Collective on Mat

6937:    Input Parameters:
6938: +  mat - the matrix
6939: .  perm - row and column permutation
6940: -  info - structure containing
6941: $      levels - number of levels of fill.
6942: $      expected fill - as ratio of original fill.

6944:    Output Parameter:
6945: .  fact - the factored matrix

6947:    Notes:
6948:    Most users should employ the KSP interface for linear solvers
6949:    instead of working directly with matrix algebra routines such as this.
6950:    See, e.g., KSPCreate().

6952:    Level: developer

6954: .seealso: MatCholeskyFactorNumeric(), MatCholeskyFactor(), MatFactorInfo

6956:     Note: this uses the definition of level of fill as in Y. Saad, 2003

6958:     Developer Note: fortran interface is not autogenerated as the f90
6959:     interface defintion cannot be generated correctly [due to MatFactorInfo]

6961:    References:
6962:      Y. Saad, Iterative methods for sparse linear systems Philadelphia: Society for Industrial and Applied Mathematics, 2003
6963: @*/
6964: PetscErrorCode MatICCFactorSymbolic(Mat fact,Mat mat,IS perm,const MatFactorInfo *info)
6965: {

6974:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
6975:   if (info->levels < 0) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_OUTOFRANGE,"Levels negative %D",(PetscInt) info->levels);
6976:   if (info->fill < 1.0) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_OUTOFRANGE,"Expected fill less than 1.0 %g",(double)info->fill);
6977:   if (!(fact)->ops->iccfactorsymbolic) {
6978:     MatSolverType stype;
6979:     MatFactorGetSolverType(fact,&stype);
6980:     SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Matrix type %s symbolic ICC using solver type %s",((PetscObject)mat)->type_name,stype);
6981:   }
6982:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
6983:   MatCheckPreallocated(mat,2);

6985:   PetscLogEventBegin(MAT_ICCFactorSymbolic,mat,perm,0,0);
6986:   (fact->ops->iccfactorsymbolic)(fact,mat,perm,info);
6987:   PetscLogEventEnd(MAT_ICCFactorSymbolic,mat,perm,0,0);
6988:   return(0);
6989: }

6991: /*@C
6992:    MatCreateSubMatrices - Extracts several submatrices from a matrix. If submat
6993:    points to an array of valid matrices, they may be reused to store the new
6994:    submatrices.

6996:    Collective on Mat

6998:    Input Parameters:
6999: +  mat - the matrix
7000: .  n   - the number of submatrixes to be extracted (on this processor, may be zero)
7001: .  irow, icol - index sets of rows and columns to extract
7002: -  scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX

7004:    Output Parameter:
7005: .  submat - the array of submatrices

7007:    Notes:
7008:    MatCreateSubMatrices() can extract ONLY sequential submatrices
7009:    (from both sequential and parallel matrices). Use MatCreateSubMatrix()
7010:    to extract a parallel submatrix.

7012:    Some matrix types place restrictions on the row and column
7013:    indices, such as that they be sorted or that they be equal to each other.

7015:    The index sets may not have duplicate entries.

7017:    When extracting submatrices from a parallel matrix, each processor can
7018:    form a different submatrix by setting the rows and columns of its
7019:    individual index sets according to the local submatrix desired.

7021:    When finished using the submatrices, the user should destroy
7022:    them with MatDestroySubMatrices().

7024:    MAT_REUSE_MATRIX can only be used when the nonzero structure of the
7025:    original matrix has not changed from that last call to MatCreateSubMatrices().

7027:    This routine creates the matrices in submat; you should NOT create them before
7028:    calling it. It also allocates the array of matrix pointers submat.

7030:    For BAIJ matrices the index sets must respect the block structure, that is if they
7031:    request one row/column in a block, they must request all rows/columns that are in
7032:    that block. For example, if the block size is 2 you cannot request just row 0 and
7033:    column 0.

7035:    Fortran Note:
7036:    The Fortran interface is slightly different from that given below; it
7037:    requires one to pass in  as submat a Mat (integer) array of size at least n+1.

7039:    Level: advanced


7042: .seealso: MatDestroySubMatrices(), MatCreateSubMatrix(), MatGetRow(), MatGetDiagonal(), MatReuse
7043: @*/
7044: PetscErrorCode MatCreateSubMatrices(Mat mat,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *submat[])
7045: {
7047:   PetscInt       i;
7048:   PetscBool      eq;

7053:   if (n) {
7058:   }
7060:   if (n && scall == MAT_REUSE_MATRIX) {
7063:   }
7064:   if (!mat->ops->createsubmatrices) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
7065:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
7066:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
7067:   MatCheckPreallocated(mat,1);

7069:   PetscLogEventBegin(MAT_CreateSubMats,mat,0,0,0);
7070:   (*mat->ops->createsubmatrices)(mat,n,irow,icol,scall,submat);
7071:   PetscLogEventEnd(MAT_CreateSubMats,mat,0,0,0);
7072:   for (i=0; i<n; i++) {
7073:     (*submat)[i]->factortype = MAT_FACTOR_NONE;  /* in case in place factorization was previously done on submatrix */
7074:     ISEqualUnsorted(irow[i],icol[i],&eq);
7075:     if (eq) {
7076:       MatPropagateSymmetryOptions(mat,(*submat)[i]);
7077:     }
7078:   }
7079:   return(0);
7080: }

7082: /*@C
7083:    MatCreateSubMatricesMPI - Extracts MPI submatrices across a sub communicator of mat (by pairs of IS that may live on subcomms).

7085:    Collective on Mat

7087:    Input Parameters:
7088: +  mat - the matrix
7089: .  n   - the number of submatrixes to be extracted
7090: .  irow, icol - index sets of rows and columns to extract
7091: -  scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX

7093:    Output Parameter:
7094: .  submat - the array of submatrices

7096:    Level: advanced


7099: .seealso: MatCreateSubMatrices(), MatCreateSubMatrix(), MatGetRow(), MatGetDiagonal(), MatReuse
7100: @*/
7101: PetscErrorCode MatCreateSubMatricesMPI(Mat mat,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *submat[])
7102: {
7104:   PetscInt       i;
7105:   PetscBool      eq;

7110:   if (n) {
7115:   }
7117:   if (n && scall == MAT_REUSE_MATRIX) {
7120:   }
7121:   if (!mat->ops->createsubmatricesmpi) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
7122:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
7123:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
7124:   MatCheckPreallocated(mat,1);

7126:   PetscLogEventBegin(MAT_CreateSubMats,mat,0,0,0);
7127:   (*mat->ops->createsubmatricesmpi)(mat,n,irow,icol,scall,submat);
7128:   PetscLogEventEnd(MAT_CreateSubMats,mat,0,0,0);
7129:   for (i=0; i<n; i++) {
7130:     ISEqualUnsorted(irow[i],icol[i],&eq);
7131:     if (eq) {
7132:       MatPropagateSymmetryOptions(mat,(*submat)[i]);
7133:     }
7134:   }
7135:   return(0);
7136: }

7138: /*@C
7139:    MatDestroyMatrices - Destroys an array of matrices.

7141:    Collective on Mat

7143:    Input Parameters:
7144: +  n - the number of local matrices
7145: -  mat - the matrices (note that this is a pointer to the array of matrices)

7147:    Level: advanced

7149:     Notes:
7150:     Frees not only the matrices, but also the array that contains the matrices
7151:            In Fortran will not free the array.

7153: .seealso: MatCreateSubMatrices() MatDestroySubMatrices()
7154: @*/
7155: PetscErrorCode MatDestroyMatrices(PetscInt n,Mat *mat[])
7156: {
7158:   PetscInt       i;

7161:   if (!*mat) return(0);
7162:   if (n < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Trying to destroy negative number of matrices %D",n);

7165:   for (i=0; i<n; i++) {
7166:     MatDestroy(&(*mat)[i]);
7167:   }

7169:   /* memory is allocated even if n = 0 */
7170:   PetscFree(*mat);
7171:   return(0);
7172: }

7174: /*@C
7175:    MatDestroySubMatrices - Destroys a set of matrices obtained with MatCreateSubMatrices().

7177:    Collective on Mat

7179:    Input Parameters:
7180: +  n - the number of local matrices
7181: -  mat - the matrices (note that this is a pointer to the array of matrices, just to match the calling
7182:                        sequence of MatCreateSubMatrices())

7184:    Level: advanced

7186:     Notes:
7187:     Frees not only the matrices, but also the array that contains the matrices
7188:            In Fortran will not free the array.

7190: .seealso: MatCreateSubMatrices()
7191: @*/
7192: PetscErrorCode MatDestroySubMatrices(PetscInt n,Mat *mat[])
7193: {
7195:   Mat            mat0;

7198:   if (!*mat) return(0);
7199:   /* mat[] is an array of length n+1, see MatCreateSubMatrices_xxx() */
7200:   if (n < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Trying to destroy negative number of matrices %D",n);

7203:   mat0 = (*mat)[0];
7204:   if (mat0 && mat0->ops->destroysubmatrices) {
7205:     (mat0->ops->destroysubmatrices)(n,mat);
7206:   } else {
7207:     MatDestroyMatrices(n,mat);
7208:   }
7209:   return(0);
7210: }

7212: /*@C
7213:    MatGetSeqNonzeroStructure - Extracts the sequential nonzero structure from a matrix.

7215:    Collective on Mat

7217:    Input Parameters:
7218: .  mat - the matrix

7220:    Output Parameter:
7221: .  matstruct - the sequential matrix with the nonzero structure of mat

7223:   Level: intermediate

7225: .seealso: MatDestroySeqNonzeroStructure(), MatCreateSubMatrices(), MatDestroyMatrices()
7226: @*/
7227: PetscErrorCode MatGetSeqNonzeroStructure(Mat mat,Mat *matstruct)
7228: {


7236:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
7237:   MatCheckPreallocated(mat,1);

7239:   if (!mat->ops->getseqnonzerostructure) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Not for matrix type %s\n",((PetscObject)mat)->type_name);
7240:   PetscLogEventBegin(MAT_GetSeqNonzeroStructure,mat,0,0,0);
7241:   (*mat->ops->getseqnonzerostructure)(mat,matstruct);
7242:   PetscLogEventEnd(MAT_GetSeqNonzeroStructure,mat,0,0,0);
7243:   return(0);
7244: }

7246: /*@C
7247:    MatDestroySeqNonzeroStructure - Destroys matrix obtained with MatGetSeqNonzeroStructure().

7249:    Collective on Mat

7251:    Input Parameters:
7252: .  mat - the matrix (note that this is a pointer to the array of matrices, just to match the calling
7253:                        sequence of MatGetSequentialNonzeroStructure())

7255:    Level: advanced

7257:     Notes:
7258:     Frees not only the matrices, but also the array that contains the matrices

7260: .seealso: MatGetSeqNonzeroStructure()
7261: @*/
7262: PetscErrorCode MatDestroySeqNonzeroStructure(Mat *mat)
7263: {

7268:   MatDestroy(mat);
7269:   return(0);
7270: }

7272: /*@
7273:    MatIncreaseOverlap - Given a set of submatrices indicated by index sets,
7274:    replaces the index sets by larger ones that represent submatrices with
7275:    additional overlap.

7277:    Collective on Mat

7279:    Input Parameters:
7280: +  mat - the matrix
7281: .  n   - the number of index sets
7282: .  is  - the array of index sets (these index sets will changed during the call)
7283: -  ov  - the additional overlap requested

7285:    Options Database:
7286: .  -mat_increase_overlap_scalable - use a scalable algorithm to compute the overlap (supported by MPIAIJ matrix)

7288:    Level: developer


7291: .seealso: MatCreateSubMatrices()
7292: @*/
7293: PetscErrorCode MatIncreaseOverlap(Mat mat,PetscInt n,IS is[],PetscInt ov)
7294: {

7300:   if (n < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must have one or more domains, you have %D",n);
7301:   if (n) {
7304:   }
7305:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
7306:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
7307:   MatCheckPreallocated(mat,1);

7309:   if (!ov) return(0);
7310:   if (!mat->ops->increaseoverlap) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
7311:   PetscLogEventBegin(MAT_IncreaseOverlap,mat,0,0,0);
7312:   (*mat->ops->increaseoverlap)(mat,n,is,ov);
7313:   PetscLogEventEnd(MAT_IncreaseOverlap,mat,0,0,0);
7314:   return(0);
7315: }


7318: PetscErrorCode MatIncreaseOverlapSplit_Single(Mat,IS*,PetscInt);

7320: /*@
7321:    MatIncreaseOverlapSplit - Given a set of submatrices indicated by index sets across
7322:    a sub communicator, replaces the index sets by larger ones that represent submatrices with
7323:    additional overlap.

7325:    Collective on Mat

7327:    Input Parameters:
7328: +  mat - the matrix
7329: .  n   - the number of index sets
7330: .  is  - the array of index sets (these index sets will changed during the call)
7331: -  ov  - the additional overlap requested

7333:    Options Database:
7334: .  -mat_increase_overlap_scalable - use a scalable algorithm to compute the overlap (supported by MPIAIJ matrix)

7336:    Level: developer


7339: .seealso: MatCreateSubMatrices()
7340: @*/
7341: PetscErrorCode MatIncreaseOverlapSplit(Mat mat,PetscInt n,IS is[],PetscInt ov)
7342: {
7343:   PetscInt       i;

7349:   if (n < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must have one or more domains, you have %D",n);
7350:   if (n) {
7353:   }
7354:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
7355:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
7356:   MatCheckPreallocated(mat,1);
7357:   if (!ov) return(0);
7358:   PetscLogEventBegin(MAT_IncreaseOverlap,mat,0,0,0);
7359:   for (i=0; i<n; i++){
7360:          MatIncreaseOverlapSplit_Single(mat,&is[i],ov);
7361:   }
7362:   PetscLogEventEnd(MAT_IncreaseOverlap,mat,0,0,0);
7363:   return(0);
7364: }




7369: /*@
7370:    MatGetBlockSize - Returns the matrix block size.

7372:    Not Collective

7374:    Input Parameter:
7375: .  mat - the matrix

7377:    Output Parameter:
7378: .  bs - block size

7380:    Notes:
7381:     Block row formats are MATSEQBAIJ, MATMPIBAIJ, MATSEQSBAIJ, MATMPISBAIJ. These formats ALWAYS have square block storage in the matrix.

7383:    If the block size has not been set yet this routine returns 1.

7385:    Level: intermediate

7387: .seealso: MatCreateSeqBAIJ(), MatCreateBAIJ(), MatGetBlockSizes()
7388: @*/
7389: PetscErrorCode MatGetBlockSize(Mat mat,PetscInt *bs)
7390: {
7394:   *bs = PetscAbs(mat->rmap->bs);
7395:   return(0);
7396: }

7398: /*@
7399:    MatGetBlockSizes - Returns the matrix block row and column sizes.

7401:    Not Collective

7403:    Input Parameter:
7404: .  mat - the matrix

7406:    Output Parameter:
7407: +  rbs - row block size
7408: -  cbs - column block size

7410:    Notes:
7411:     Block row formats are MATSEQBAIJ, MATMPIBAIJ, MATSEQSBAIJ, MATMPISBAIJ. These formats ALWAYS have square block storage in the matrix.
7412:     If you pass a different block size for the columns than the rows, the row block size determines the square block storage.

7414:    If a block size has not been set yet this routine returns 1.

7416:    Level: intermediate

7418: .seealso: MatCreateSeqBAIJ(), MatCreateBAIJ(), MatGetBlockSize(), MatSetBlockSize(), MatSetBlockSizes()
7419: @*/
7420: PetscErrorCode MatGetBlockSizes(Mat mat,PetscInt *rbs, PetscInt *cbs)
7421: {
7426:   if (rbs) *rbs = PetscAbs(mat->rmap->bs);
7427:   if (cbs) *cbs = PetscAbs(mat->cmap->bs);
7428:   return(0);
7429: }

7431: /*@
7432:    MatSetBlockSize - Sets the matrix block size.

7434:    Logically Collective on Mat

7436:    Input Parameters:
7437: +  mat - the matrix
7438: -  bs - block size

7440:    Notes:
7441:     Block row formats are MATSEQBAIJ, MATMPIBAIJ, MATSEQSBAIJ, MATMPISBAIJ. These formats ALWAYS have square block storage in the matrix.
7442:     This must be called before MatSetUp() or MatXXXSetPreallocation() (or will default to 1) and the block size cannot be changed later.

7444:     For MATMPIAIJ and MATSEQAIJ matrix formats, this function can be called at a later stage, provided that the specified block size
7445:     is compatible with the matrix local sizes.

7447:    Level: intermediate

7449: .seealso: MatCreateSeqBAIJ(), MatCreateBAIJ(), MatGetBlockSize(), MatSetBlockSizes(), MatGetBlockSizes()
7450: @*/
7451: PetscErrorCode MatSetBlockSize(Mat mat,PetscInt bs)
7452: {

7458:   MatSetBlockSizes(mat,bs,bs);
7459:   return(0);
7460: }

7462: /*@
7463:    MatSetVariableBlockSizes - Sets a diagonal blocks of the matrix that need not be of the same size

7465:    Logically Collective on Mat

7467:    Input Parameters:
7468: +  mat - the matrix
7469: .  nblocks - the number of blocks on this process
7470: -  bsizes - the block sizes

7472:    Notes:
7473:     Currently used by PCVPBJACOBI for SeqAIJ matrices

7475:    Level: intermediate

7477: .seealso: MatCreateSeqBAIJ(), MatCreateBAIJ(), MatGetBlockSize(), MatSetBlockSizes(), MatGetBlockSizes(), MatGetVariableBlockSizes()
7478: @*/
7479: PetscErrorCode MatSetVariableBlockSizes(Mat mat,PetscInt nblocks,PetscInt *bsizes)
7480: {
7482:   PetscInt       i,ncnt = 0, nlocal;

7486:   if (nblocks < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number of local blocks must be great than or equal to zero");
7487:   MatGetLocalSize(mat,&nlocal,NULL);
7488:   for (i=0; i<nblocks; i++) ncnt += bsizes[i];
7489:   if (ncnt != nlocal) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Sum of local block sizes %D does not equal local size of matrix %D",ncnt,nlocal);
7490:   PetscFree(mat->bsizes);
7491:   mat->nblocks = nblocks;
7492:   PetscMalloc1(nblocks,&mat->bsizes);
7493:   PetscArraycpy(mat->bsizes,bsizes,nblocks);
7494:   return(0);
7495: }

7497: /*@C
7498:    MatGetVariableBlockSizes - Gets a diagonal blocks of the matrix that need not be of the same size

7500:    Logically Collective on Mat

7502:    Input Parameters:
7503: .  mat - the matrix

7505:    Output Parameters:
7506: +  nblocks - the number of blocks on this process
7507: -  bsizes - the block sizes

7509:    Notes: Currently not supported from Fortran

7511:    Level: intermediate

7513: .seealso: MatCreateSeqBAIJ(), MatCreateBAIJ(), MatGetBlockSize(), MatSetBlockSizes(), MatGetBlockSizes(), MatSetVariableBlockSizes()
7514: @*/
7515: PetscErrorCode MatGetVariableBlockSizes(Mat mat,PetscInt *nblocks,const PetscInt **bsizes)
7516: {
7519:   *nblocks = mat->nblocks;
7520:   *bsizes  = mat->bsizes;
7521:   return(0);
7522: }

7524: /*@
7525:    MatSetBlockSizes - Sets the matrix block row and column sizes.

7527:    Logically Collective on Mat

7529:    Input Parameters:
7530: +  mat - the matrix
7531: .  rbs - row block size
7532: -  cbs - column block size

7534:    Notes:
7535:     Block row formats are MATSEQBAIJ, MATMPIBAIJ, MATSEQSBAIJ, MATMPISBAIJ. These formats ALWAYS have square block storage in the matrix.
7536:     If you pass a different block size for the columns than the rows, the row block size determines the square block storage.
7537:     This must be called before MatSetUp() or MatXXXSetPreallocation() (or will default to 1) and the block size cannot be changed later.

7539:     For MATMPIAIJ and MATSEQAIJ matrix formats, this function can be called at a later stage, provided that the specified block sizes
7540:     are compatible with the matrix local sizes.

7542:     The row and column block size determine the blocksize of the "row" and "column" vectors returned by MatCreateVecs().

7544:    Level: intermediate

7546: .seealso: MatCreateSeqBAIJ(), MatCreateBAIJ(), MatGetBlockSize(), MatSetBlockSize(), MatGetBlockSizes()
7547: @*/
7548: PetscErrorCode MatSetBlockSizes(Mat mat,PetscInt rbs,PetscInt cbs)
7549: {

7556:   if (mat->ops->setblocksizes) {
7557:     (*mat->ops->setblocksizes)(mat,rbs,cbs);
7558:   }
7559:   if (mat->rmap->refcnt) {
7560:     ISLocalToGlobalMapping l2g = NULL;
7561:     PetscLayout            nmap = NULL;

7563:     PetscLayoutDuplicate(mat->rmap,&nmap);
7564:     if (mat->rmap->mapping) {
7565:       ISLocalToGlobalMappingDuplicate(mat->rmap->mapping,&l2g);
7566:     }
7567:     PetscLayoutDestroy(&mat->rmap);
7568:     mat->rmap = nmap;
7569:     mat->rmap->mapping = l2g;
7570:   }
7571:   if (mat->cmap->refcnt) {
7572:     ISLocalToGlobalMapping l2g = NULL;
7573:     PetscLayout            nmap = NULL;

7575:     PetscLayoutDuplicate(mat->cmap,&nmap);
7576:     if (mat->cmap->mapping) {
7577:       ISLocalToGlobalMappingDuplicate(mat->cmap->mapping,&l2g);
7578:     }
7579:     PetscLayoutDestroy(&mat->cmap);
7580:     mat->cmap = nmap;
7581:     mat->cmap->mapping = l2g;
7582:   }
7583:   PetscLayoutSetBlockSize(mat->rmap,rbs);
7584:   PetscLayoutSetBlockSize(mat->cmap,cbs);
7585:   return(0);
7586: }

7588: /*@
7589:    MatSetBlockSizesFromMats - Sets the matrix block row and column sizes to match a pair of matrices

7591:    Logically Collective on Mat

7593:    Input Parameters:
7594: +  mat - the matrix
7595: .  fromRow - matrix from which to copy row block size
7596: -  fromCol - matrix from which to copy column block size (can be same as fromRow)

7598:    Level: developer

7600: .seealso: MatCreateSeqBAIJ(), MatCreateBAIJ(), MatGetBlockSize(), MatSetBlockSizes()
7601: @*/
7602: PetscErrorCode MatSetBlockSizesFromMats(Mat mat,Mat fromRow,Mat fromCol)
7603: {

7610:   if (fromRow->rmap->bs > 0) {PetscLayoutSetBlockSize(mat->rmap,fromRow->rmap->bs);}
7611:   if (fromCol->cmap->bs > 0) {PetscLayoutSetBlockSize(mat->cmap,fromCol->cmap->bs);}
7612:   return(0);
7613: }

7615: /*@
7616:    MatResidual - Default routine to calculate the residual.

7618:    Collective on Mat

7620:    Input Parameters:
7621: +  mat - the matrix
7622: .  b   - the right-hand-side
7623: -  x   - the approximate solution

7625:    Output Parameter:
7626: .  r - location to store the residual

7628:    Level: developer

7630: .seealso: PCMGSetResidual()
7631: @*/
7632: PetscErrorCode MatResidual(Mat mat,Vec b,Vec x,Vec r)
7633: {

7642:   MatCheckPreallocated(mat,1);
7643:   PetscLogEventBegin(MAT_Residual,mat,0,0,0);
7644:   if (!mat->ops->residual) {
7645:     MatMult(mat,x,r);
7646:     VecAYPX(r,-1.0,b);
7647:   } else {
7648:     (*mat->ops->residual)(mat,b,x,r);
7649:   }
7650:   PetscLogEventEnd(MAT_Residual,mat,0,0,0);
7651:   return(0);
7652: }

7654: /*@C
7655:     MatGetRowIJ - Returns the compressed row storage i and j indices for sequential matrices.

7657:    Collective on Mat

7659:     Input Parameters:
7660: +   mat - the matrix
7661: .   shift -  0 or 1 indicating we want the indices starting at 0 or 1
7662: .   symmetric - PETSC_TRUE or PETSC_FALSE indicating the matrix data structure should be   symmetrized
7663: -   inodecompressed - PETSC_TRUE or PETSC_FALSE  indicating if the nonzero structure of the
7664:                  inodes or the nonzero elements is wanted. For BAIJ matrices the compressed version is
7665:                  always used.

7667:     Output Parameters:
7668: +   n - number of rows in the (possibly compressed) matrix
7669: .   ia - the row pointers; that is ia[0] = 0, ia[row] = ia[row-1] + number of elements in that row of the matrix
7670: .   ja - the column indices
7671: -   done - indicates if the routine actually worked and returned appropriate ia[] and ja[] arrays; callers
7672:            are responsible for handling the case when done == PETSC_FALSE and ia and ja are not set

7674:     Level: developer

7676:     Notes:
7677:     You CANNOT change any of the ia[] or ja[] values.

7679:     Use MatRestoreRowIJ() when you are finished accessing the ia[] and ja[] values.

7681:     Fortran Notes:
7682:     In Fortran use
7683: $
7684: $      PetscInt ia(1), ja(1)
7685: $      PetscOffset iia, jja
7686: $      call MatGetRowIJ(mat,shift,symmetric,inodecompressed,n,ia,iia,ja,jja,done,ierr)
7687: $      ! Access the ith and jth entries via ia(iia + i) and ja(jja + j)

7689:      or
7690: $
7691: $    PetscInt, pointer :: ia(:),ja(:)
7692: $    call MatGetRowIJF90(mat,shift,symmetric,inodecompressed,n,ia,ja,done,ierr)
7693: $    ! Access the ith and jth entries via ia(i) and ja(j)

7695: .seealso: MatGetColumnIJ(), MatRestoreRowIJ(), MatSeqAIJGetArray()
7696: @*/
7697: PetscErrorCode MatGetRowIJ(Mat mat,PetscInt shift,PetscBool symmetric,PetscBool inodecompressed,PetscInt *n,const PetscInt *ia[],const PetscInt *ja[],PetscBool  *done)
7698: {

7708:   MatCheckPreallocated(mat,1);
7709:   if (!mat->ops->getrowij) *done = PETSC_FALSE;
7710:   else {
7711:     *done = PETSC_TRUE;
7712:     PetscLogEventBegin(MAT_GetRowIJ,mat,0,0,0);
7713:     (*mat->ops->getrowij)(mat,shift,symmetric,inodecompressed,n,ia,ja,done);
7714:     PetscLogEventEnd(MAT_GetRowIJ,mat,0,0,0);
7715:   }
7716:   return(0);
7717: }

7719: /*@C
7720:     MatGetColumnIJ - Returns the compressed column storage i and j indices for sequential matrices.

7722:     Collective on Mat

7724:     Input Parameters:
7725: +   mat - the matrix
7726: .   shift - 1 or zero indicating we want the indices starting at 0 or 1
7727: .   symmetric - PETSC_TRUE or PETSC_FALSE indicating the matrix data structure should be
7728:                 symmetrized
7729: .   inodecompressed - PETSC_TRUE or PETSC_FALSE indicating if the nonzero structure of the
7730:                  inodes or the nonzero elements is wanted. For BAIJ matrices the compressed version is
7731:                  always used.
7732: .   n - number of columns in the (possibly compressed) matrix
7733: .   ia - the column pointers; that is ia[0] = 0, ia[col] = i[col-1] + number of elements in that col of the matrix
7734: -   ja - the row indices

7736:     Output Parameters:
7737: .   done - PETSC_TRUE or PETSC_FALSE, indicating whether the values have been returned

7739:     Level: developer

7741: .seealso: MatGetRowIJ(), MatRestoreColumnIJ()
7742: @*/
7743: PetscErrorCode MatGetColumnIJ(Mat mat,PetscInt shift,PetscBool symmetric,PetscBool inodecompressed,PetscInt *n,const PetscInt *ia[],const PetscInt *ja[],PetscBool  *done)
7744: {

7754:   MatCheckPreallocated(mat,1);
7755:   if (!mat->ops->getcolumnij) *done = PETSC_FALSE;
7756:   else {
7757:     *done = PETSC_TRUE;
7758:     (*mat->ops->getcolumnij)(mat,shift,symmetric,inodecompressed,n,ia,ja,done);
7759:   }
7760:   return(0);
7761: }

7763: /*@C
7764:     MatRestoreRowIJ - Call after you are completed with the ia,ja indices obtained with
7765:     MatGetRowIJ().

7767:     Collective on Mat

7769:     Input Parameters:
7770: +   mat - the matrix
7771: .   shift - 1 or zero indicating we want the indices starting at 0 or 1
7772: .   symmetric - PETSC_TRUE or PETSC_FALSE indicating the matrix data structure should be
7773:                 symmetrized
7774: .   inodecompressed -  PETSC_TRUE or PETSC_FALSE indicating if the nonzero structure of the
7775:                  inodes or the nonzero elements is wanted. For BAIJ matrices the compressed version is
7776:                  always used.
7777: .   n - size of (possibly compressed) matrix
7778: .   ia - the row pointers
7779: -   ja - the column indices

7781:     Output Parameters:
7782: .   done - PETSC_TRUE or PETSC_FALSE indicated that the values have been returned

7784:     Note:
7785:     This routine zeros out n, ia, and ja. This is to prevent accidental
7786:     us of the array after it has been restored. If you pass NULL, it will
7787:     not zero the pointers.  Use of ia or ja after MatRestoreRowIJ() is invalid.

7789:     Level: developer

7791: .seealso: MatGetRowIJ(), MatRestoreColumnIJ()
7792: @*/
7793: PetscErrorCode MatRestoreRowIJ(Mat mat,PetscInt shift,PetscBool symmetric,PetscBool inodecompressed,PetscInt *n,const PetscInt *ia[],const PetscInt *ja[],PetscBool  *done)
7794: {

7803:   MatCheckPreallocated(mat,1);

7805:   if (!mat->ops->restorerowij) *done = PETSC_FALSE;
7806:   else {
7807:     *done = PETSC_TRUE;
7808:     (*mat->ops->restorerowij)(mat,shift,symmetric,inodecompressed,n,ia,ja,done);
7809:     if (n)  *n = 0;
7810:     if (ia) *ia = NULL;
7811:     if (ja) *ja = NULL;
7812:   }
7813:   return(0);
7814: }

7816: /*@C
7817:     MatRestoreColumnIJ - Call after you are completed with the ia,ja indices obtained with
7818:     MatGetColumnIJ().

7820:     Collective on Mat

7822:     Input Parameters:
7823: +   mat - the matrix
7824: .   shift - 1 or zero indicating we want the indices starting at 0 or 1
7825: -   symmetric - PETSC_TRUE or PETSC_FALSE indicating the matrix data structure should be
7826:                 symmetrized
7827: -   inodecompressed - PETSC_TRUE or PETSC_FALSE indicating if the nonzero structure of the
7828:                  inodes or the nonzero elements is wanted. For BAIJ matrices the compressed version is
7829:                  always used.

7831:     Output Parameters:
7832: +   n - size of (possibly compressed) matrix
7833: .   ia - the column pointers
7834: .   ja - the row indices
7835: -   done - PETSC_TRUE or PETSC_FALSE indicated that the values have been returned

7837:     Level: developer

7839: .seealso: MatGetColumnIJ(), MatRestoreRowIJ()
7840: @*/
7841: PetscErrorCode MatRestoreColumnIJ(Mat mat,PetscInt shift,PetscBool symmetric,PetscBool inodecompressed,PetscInt *n,const PetscInt *ia[],const PetscInt *ja[],PetscBool  *done)
7842: {

7851:   MatCheckPreallocated(mat,1);

7853:   if (!mat->ops->restorecolumnij) *done = PETSC_FALSE;
7854:   else {
7855:     *done = PETSC_TRUE;
7856:     (*mat->ops->restorecolumnij)(mat,shift,symmetric,inodecompressed,n,ia,ja,done);
7857:     if (n)  *n = 0;
7858:     if (ia) *ia = NULL;
7859:     if (ja) *ja = NULL;
7860:   }
7861:   return(0);
7862: }

7864: /*@C
7865:     MatColoringPatch -Used inside matrix coloring routines that
7866:     use MatGetRowIJ() and/or MatGetColumnIJ().

7868:     Collective on Mat

7870:     Input Parameters:
7871: +   mat - the matrix
7872: .   ncolors - max color value
7873: .   n   - number of entries in colorarray
7874: -   colorarray - array indicating color for each column

7876:     Output Parameters:
7877: .   iscoloring - coloring generated using colorarray information

7879:     Level: developer

7881: .seealso: MatGetRowIJ(), MatGetColumnIJ()

7883: @*/
7884: PetscErrorCode MatColoringPatch(Mat mat,PetscInt ncolors,PetscInt n,ISColoringValue colorarray[],ISColoring *iscoloring)
7885: {

7893:   MatCheckPreallocated(mat,1);

7895:   if (!mat->ops->coloringpatch) {
7896:     ISColoringCreate(PetscObjectComm((PetscObject)mat),ncolors,n,colorarray,PETSC_OWN_POINTER,iscoloring);
7897:   } else {
7898:     (*mat->ops->coloringpatch)(mat,ncolors,n,colorarray,iscoloring);
7899:   }
7900:   return(0);
7901: }


7904: /*@
7905:    MatSetUnfactored - Resets a factored matrix to be treated as unfactored.

7907:    Logically Collective on Mat

7909:    Input Parameter:
7910: .  mat - the factored matrix to be reset

7912:    Notes:
7913:    This routine should be used only with factored matrices formed by in-place
7914:    factorization via ILU(0) (or by in-place LU factorization for the MATSEQDENSE
7915:    format).  This option can save memory, for example, when solving nonlinear
7916:    systems with a matrix-free Newton-Krylov method and a matrix-based, in-place
7917:    ILU(0) preconditioner.

7919:    Note that one can specify in-place ILU(0) factorization by calling
7920: .vb
7921:      PCType(pc,PCILU);
7922:      PCFactorSeUseInPlace(pc);
7923: .ve
7924:    or by using the options -pc_type ilu -pc_factor_in_place

7926:    In-place factorization ILU(0) can also be used as a local
7927:    solver for the blocks within the block Jacobi or additive Schwarz
7928:    methods (runtime option: -sub_pc_factor_in_place).  See Users-Manual: ch_pc
7929:    for details on setting local solver options.

7931:    Most users should employ the simplified KSP interface for linear solvers
7932:    instead of working directly with matrix algebra routines such as this.
7933:    See, e.g., KSPCreate().

7935:    Level: developer

7937: .seealso: PCFactorSetUseInPlace(), PCFactorGetUseInPlace()

7939: @*/
7940: PetscErrorCode MatSetUnfactored(Mat mat)
7941: {

7947:   MatCheckPreallocated(mat,1);
7948:   mat->factortype = MAT_FACTOR_NONE;
7949:   if (!mat->ops->setunfactored) return(0);
7950:   (*mat->ops->setunfactored)(mat);
7951:   return(0);
7952: }

7954: /*MC
7955:     MatDenseGetArrayF90 - Accesses a matrix array from Fortran90.

7957:     Synopsis:
7958:     MatDenseGetArrayF90(Mat x,{Scalar, pointer :: xx_v(:,:)},integer ierr)

7960:     Not collective

7962:     Input Parameter:
7963: .   x - matrix

7965:     Output Parameters:
7966: +   xx_v - the Fortran90 pointer to the array
7967: -   ierr - error code

7969:     Example of Usage:
7970: .vb
7971:       PetscScalar, pointer xx_v(:,:)
7972:       ....
7973:       call MatDenseGetArrayF90(x,xx_v,ierr)
7974:       a = xx_v(3)
7975:       call MatDenseRestoreArrayF90(x,xx_v,ierr)
7976: .ve

7978:     Level: advanced

7980: .seealso:  MatDenseRestoreArrayF90(), MatDenseGetArray(), MatDenseRestoreArray(), MatSeqAIJGetArrayF90()

7982: M*/

7984: /*MC
7985:     MatDenseRestoreArrayF90 - Restores a matrix array that has been
7986:     accessed with MatDenseGetArrayF90().

7988:     Synopsis:
7989:     MatDenseRestoreArrayF90(Mat x,{Scalar, pointer :: xx_v(:,:)},integer ierr)

7991:     Not collective

7993:     Input Parameters:
7994: +   x - matrix
7995: -   xx_v - the Fortran90 pointer to the array

7997:     Output Parameter:
7998: .   ierr - error code

8000:     Example of Usage:
8001: .vb
8002:        PetscScalar, pointer xx_v(:,:)
8003:        ....
8004:        call MatDenseGetArrayF90(x,xx_v,ierr)
8005:        a = xx_v(3)
8006:        call MatDenseRestoreArrayF90(x,xx_v,ierr)
8007: .ve

8009:     Level: advanced

8011: .seealso:  MatDenseGetArrayF90(), MatDenseGetArray(), MatDenseRestoreArray(), MatSeqAIJRestoreArrayF90()

8013: M*/


8016: /*MC
8017:     MatSeqAIJGetArrayF90 - Accesses a matrix array from Fortran90.

8019:     Synopsis:
8020:     MatSeqAIJGetArrayF90(Mat x,{Scalar, pointer :: xx_v(:)},integer ierr)

8022:     Not collective

8024:     Input Parameter:
8025: .   x - matrix

8027:     Output Parameters:
8028: +   xx_v - the Fortran90 pointer to the array
8029: -   ierr - error code

8031:     Example of Usage:
8032: .vb
8033:       PetscScalar, pointer xx_v(:)
8034:       ....
8035:       call MatSeqAIJGetArrayF90(x,xx_v,ierr)
8036:       a = xx_v(3)
8037:       call MatSeqAIJRestoreArrayF90(x,xx_v,ierr)
8038: .ve

8040:     Level: advanced

8042: .seealso:  MatSeqAIJRestoreArrayF90(), MatSeqAIJGetArray(), MatSeqAIJRestoreArray(), MatDenseGetArrayF90()

8044: M*/

8046: /*MC
8047:     MatSeqAIJRestoreArrayF90 - Restores a matrix array that has been
8048:     accessed with MatSeqAIJGetArrayF90().

8050:     Synopsis:
8051:     MatSeqAIJRestoreArrayF90(Mat x,{Scalar, pointer :: xx_v(:)},integer ierr)

8053:     Not collective

8055:     Input Parameters:
8056: +   x - matrix
8057: -   xx_v - the Fortran90 pointer to the array

8059:     Output Parameter:
8060: .   ierr - error code

8062:     Example of Usage:
8063: .vb
8064:        PetscScalar, pointer xx_v(:)
8065:        ....
8066:        call MatSeqAIJGetArrayF90(x,xx_v,ierr)
8067:        a = xx_v(3)
8068:        call MatSeqAIJRestoreArrayF90(x,xx_v,ierr)
8069: .ve

8071:     Level: advanced

8073: .seealso:  MatSeqAIJGetArrayF90(), MatSeqAIJGetArray(), MatSeqAIJRestoreArray(), MatDenseRestoreArrayF90()

8075: M*/


8078: /*@
8079:     MatCreateSubMatrix - Gets a single submatrix on the same number of processors
8080:                       as the original matrix.

8082:     Collective on Mat

8084:     Input Parameters:
8085: +   mat - the original matrix
8086: .   isrow - parallel IS containing the rows this processor should obtain
8087: .   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.
8088: -   cll - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX

8090:     Output Parameter:
8091: .   newmat - the new submatrix, of the same type as the old

8093:     Level: advanced

8095:     Notes:
8096:     The submatrix will be able to be multiplied with vectors using the same layout as iscol.

8098:     Some matrix types place restrictions on the row and column indices, such
8099:     as that they be sorted or that they be equal to each other.

8101:     The index sets may not have duplicate entries.

8103:       The first time this is called you should use a cll of MAT_INITIAL_MATRIX,
8104:    the MatCreateSubMatrix() routine will create the newmat for you. Any additional calls
8105:    to this routine with a mat of the same nonzero structure and with a call of MAT_REUSE_MATRIX
8106:    will reuse the matrix generated the first time.  You should call MatDestroy() on newmat when
8107:    you are finished using it.

8109:     The communicator of the newly obtained matrix is ALWAYS the same as the communicator of
8110:     the input matrix.

8112:     If iscol is NULL then all columns are obtained (not supported in Fortran).

8114:    Example usage:
8115:    Consider the following 8x8 matrix with 34 non-zero values, that is
8116:    assembled across 3 processors. Let's assume that proc0 owns 3 rows,
8117:    proc1 owns 3 rows, proc2 owns 2 rows. This division can be shown
8118:    as follows:

8120: .vb
8121:             1  2  0  |  0  3  0  |  0  4
8122:     Proc0   0  5  6  |  7  0  0  |  8  0
8123:             9  0 10  | 11  0  0  | 12  0
8124:     -------------------------------------
8125:            13  0 14  | 15 16 17  |  0  0
8126:     Proc1   0 18  0  | 19 20 21  |  0  0
8127:             0  0  0  | 22 23  0  | 24  0
8128:     -------------------------------------
8129:     Proc2  25 26 27  |  0  0 28  | 29  0
8130:            30  0  0  | 31 32 33  |  0 34
8131: .ve

8133:     Suppose isrow = [0 1 | 4 | 6 7] and iscol = [1 2 | 3 4 5 | 6].  The resulting submatrix is

8135: .vb
8136:             2  0  |  0  3  0  |  0
8137:     Proc0   5  6  |  7  0  0  |  8
8138:     -------------------------------
8139:     Proc1  18  0  | 19 20 21  |  0
8140:     -------------------------------
8141:     Proc2  26 27  |  0  0 28  | 29
8142:             0  0  | 31 32 33  |  0
8143: .ve


8146: .seealso: MatCreateSubMatrices(), MatCreateSubMatricesMPI(), MatCreateSubMatrixVirtual(), MatSubMatrixVirtualUpdate()
8147: @*/
8148: PetscErrorCode MatCreateSubMatrix(Mat mat,IS isrow,IS iscol,MatReuse cll,Mat *newmat)
8149: {
8151:   PetscMPIInt    size;
8152:   Mat            *local;
8153:   IS             iscoltmp;
8154:   PetscBool      flg;

8163:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
8164:   if (cll == MAT_IGNORE_MATRIX) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Cannot use MAT_IGNORE_MATRIX");

8166:   MatCheckPreallocated(mat,1);
8167:   MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);

8169:   if (!iscol || isrow == iscol) {
8170:     PetscBool   stride;
8171:     PetscMPIInt grabentirematrix = 0,grab;
8172:     PetscObjectTypeCompare((PetscObject)isrow,ISSTRIDE,&stride);
8173:     if (stride) {
8174:       PetscInt first,step,n,rstart,rend;
8175:       ISStrideGetInfo(isrow,&first,&step);
8176:       if (step == 1) {
8177:         MatGetOwnershipRange(mat,&rstart,&rend);
8178:         if (rstart == first) {
8179:           ISGetLocalSize(isrow,&n);
8180:           if (n == rend-rstart) {
8181:             grabentirematrix = 1;
8182:           }
8183:         }
8184:       }
8185:     }
8186:     MPIU_Allreduce(&grabentirematrix,&grab,1,MPI_INT,MPI_MIN,PetscObjectComm((PetscObject)mat));
8187:     if (grab) {
8188:       PetscInfo(mat,"Getting entire matrix as submatrix\n");
8189:       if (cll == MAT_INITIAL_MATRIX) {
8190:         *newmat = mat;
8191:         PetscObjectReference((PetscObject)mat);
8192:       }
8193:       return(0);
8194:     }
8195:   }

8197:   if (!iscol) {
8198:     ISCreateStride(PetscObjectComm((PetscObject)mat),mat->cmap->n,mat->cmap->rstart,1,&iscoltmp);
8199:   } else {
8200:     iscoltmp = iscol;
8201:   }

8203:   /* if original matrix is on just one processor then use submatrix generated */
8204:   if (mat->ops->createsubmatrices && !mat->ops->createsubmatrix && size == 1 && cll == MAT_REUSE_MATRIX) {
8205:     MatCreateSubMatrices(mat,1,&isrow,&iscoltmp,MAT_REUSE_MATRIX,&newmat);
8206:     goto setproperties;
8207:   } else if (mat->ops->createsubmatrices && !mat->ops->createsubmatrix && size == 1) {
8208:     MatCreateSubMatrices(mat,1,&isrow,&iscoltmp,MAT_INITIAL_MATRIX,&local);
8209:     *newmat = *local;
8210:     PetscFree(local);
8211:     goto setproperties;
8212:   } else if (!mat->ops->createsubmatrix) {
8213:     /* Create a new matrix type that implements the operation using the full matrix */
8214:     PetscLogEventBegin(MAT_CreateSubMat,mat,0,0,0);
8215:     switch (cll) {
8216:     case MAT_INITIAL_MATRIX:
8217:       MatCreateSubMatrixVirtual(mat,isrow,iscoltmp,newmat);
8218:       break;
8219:     case MAT_REUSE_MATRIX:
8220:       MatSubMatrixVirtualUpdate(*newmat,mat,isrow,iscoltmp);
8221:       break;
8222:     default: SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_OUTOFRANGE,"Invalid MatReuse, must be either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX");
8223:     }
8224:     PetscLogEventEnd(MAT_CreateSubMat,mat,0,0,0);
8225:     goto setproperties;
8226:   }

8228:   if (!mat->ops->createsubmatrix) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
8229:   PetscLogEventBegin(MAT_CreateSubMat,mat,0,0,0);
8230:   (*mat->ops->createsubmatrix)(mat,isrow,iscoltmp,cll,newmat);
8231:   PetscLogEventEnd(MAT_CreateSubMat,mat,0,0,0);

8233: setproperties:
8234:   ISEqualUnsorted(isrow,iscoltmp,&flg);
8235:   if (flg) {
8236:     MatPropagateSymmetryOptions(mat,*newmat);
8237:   }
8238:   if (!iscol) {ISDestroy(&iscoltmp);}
8239:   if (*newmat && cll == MAT_INITIAL_MATRIX) {PetscObjectStateIncrease((PetscObject)*newmat);}
8240:   return(0);
8241: }

8243: /*@
8244:    MatPropagateSymmetryOptions - Propagates symmetry options set on a matrix to another matrix

8246:    Not Collective

8248:    Input Parameters:
8249: +  A - the matrix we wish to propagate options from
8250: -  B - the matrix we wish to propagate options to

8252:    Level: beginner

8254:    Notes: Propagates the options associated to MAT_SYMMETRY_ETERNAL, MAT_STRUCTURALLY_SYMMETRIC, MAT_HERMITIAN, MAT_SPD and MAT_SYMMETRIC

8256: .seealso: MatSetOption()
8257: @*/
8258: PetscErrorCode MatPropagateSymmetryOptions(Mat A, Mat B)
8259: {

8265:   if (A->symmetric_eternal) { /* symmetric_eternal does not have a corresponding *set flag */
8266:     MatSetOption(B,MAT_SYMMETRY_ETERNAL,A->symmetric_eternal);
8267:   }
8268:   if (A->structurally_symmetric_set) {
8269:     MatSetOption(B,MAT_STRUCTURALLY_SYMMETRIC,A->structurally_symmetric);
8270:   }
8271:   if (A->hermitian_set) {
8272:     MatSetOption(B,MAT_HERMITIAN,A->hermitian);
8273:   }
8274:   if (A->spd_set) {
8275:     MatSetOption(B,MAT_SPD,A->spd);
8276:   }
8277:   if (A->symmetric_set) {
8278:     MatSetOption(B,MAT_SYMMETRIC,A->symmetric);
8279:   }
8280:   return(0);
8281: }

8283: /*@
8284:    MatStashSetInitialSize - sets the sizes of the matrix stash, that is
8285:    used during the assembly process to store values that belong to
8286:    other processors.

8288:    Not Collective

8290:    Input Parameters:
8291: +  mat   - the matrix
8292: .  size  - the initial size of the stash.
8293: -  bsize - the initial size of the block-stash(if used).

8295:    Options Database Keys:
8296: +   -matstash_initial_size <size> or <size0,size1,...sizep-1>
8297: -   -matstash_block_initial_size <bsize>  or <bsize0,bsize1,...bsizep-1>

8299:    Level: intermediate

8301:    Notes:
8302:      The block-stash is used for values set with MatSetValuesBlocked() while
8303:      the stash is used for values set with MatSetValues()

8305:      Run with the option -info and look for output of the form
8306:      MatAssemblyBegin_MPIXXX:Stash has MM entries, uses nn mallocs.
8307:      to determine the appropriate value, MM, to use for size and
8308:      MatAssemblyBegin_MPIXXX:Block-Stash has BMM entries, uses nn mallocs.
8309:      to determine the value, BMM to use for bsize


8312: .seealso: MatAssemblyBegin(), MatAssemblyEnd(), Mat, MatStashGetInfo()

8314: @*/
8315: PetscErrorCode MatStashSetInitialSize(Mat mat,PetscInt size, PetscInt bsize)
8316: {

8322:   MatStashSetInitialSize_Private(&mat->stash,size);
8323:   MatStashSetInitialSize_Private(&mat->bstash,bsize);
8324:   return(0);
8325: }

8327: /*@
8328:    MatInterpolateAdd - w = y + A*x or A'*x depending on the shape of
8329:      the matrix

8331:    Neighbor-wise Collective on Mat

8333:    Input Parameters:
8334: +  mat   - the matrix
8335: .  x,y - the vectors
8336: -  w - where the result is stored

8338:    Level: intermediate

8340:    Notes:
8341:     w may be the same vector as y.

8343:     This allows one to use either the restriction or interpolation (its transpose)
8344:     matrix to do the interpolation

8346: .seealso: MatMultAdd(), MatMultTransposeAdd(), MatRestrict()

8348: @*/
8349: PetscErrorCode MatInterpolateAdd(Mat A,Vec x,Vec y,Vec w)
8350: {
8352:   PetscInt       M,N,Ny;

8359:   MatGetSize(A,&M,&N);
8360:   VecGetSize(y,&Ny);
8361:   if (M == Ny) {
8362:     MatMultAdd(A,x,y,w);
8363:   } else {
8364:     MatMultTransposeAdd(A,x,y,w);
8365:   }
8366:   return(0);
8367: }

8369: /*@
8370:    MatInterpolate - y = A*x or A'*x depending on the shape of
8371:      the matrix

8373:    Neighbor-wise Collective on Mat

8375:    Input Parameters:
8376: +  mat   - the matrix
8377: -  x,y - the vectors

8379:    Level: intermediate

8381:    Notes:
8382:     This allows one to use either the restriction or interpolation (its transpose)
8383:     matrix to do the interpolation

8385: .seealso: MatMultAdd(), MatMultTransposeAdd(), MatRestrict()

8387: @*/
8388: PetscErrorCode MatInterpolate(Mat A,Vec x,Vec y)
8389: {
8391:   PetscInt       M,N,Ny;

8397:   MatGetSize(A,&M,&N);
8398:   VecGetSize(y,&Ny);
8399:   if (M == Ny) {
8400:     MatMult(A,x,y);
8401:   } else {
8402:     MatMultTranspose(A,x,y);
8403:   }
8404:   return(0);
8405: }

8407: /*@
8408:    MatRestrict - y = A*x or A'*x

8410:    Neighbor-wise Collective on Mat

8412:    Input Parameters:
8413: +  mat   - the matrix
8414: -  x,y - the vectors

8416:    Level: intermediate

8418:    Notes:
8419:     This allows one to use either the restriction or interpolation (its transpose)
8420:     matrix to do the restriction

8422: .seealso: MatMultAdd(), MatMultTransposeAdd(), MatInterpolate()

8424: @*/
8425: PetscErrorCode MatRestrict(Mat A,Vec x,Vec y)
8426: {
8428:   PetscInt       M,N,Ny;

8434:   MatGetSize(A,&M,&N);
8435:   VecGetSize(y,&Ny);
8436:   if (M == Ny) {
8437:     MatMult(A,x,y);
8438:   } else {
8439:     MatMultTranspose(A,x,y);
8440:   }
8441:   return(0);
8442: }

8444: /*@
8445:    MatMatInterpolateAdd - Y = W + A*X or W + A'*X

8447:    Neighbor-wise Collective on Mat

8449:    Input Parameters:
8450: +  mat   - the matrix
8451: -  w, x - the input dense matrices

8453:    Output Parameters:
8454: .  y - the output dense matrix

8456:    Level: intermediate

8458:    Notes:
8459:     This allows one to use either the restriction or interpolation (its transpose)
8460:     matrix to do the interpolation. y matrix can be reused if already created with the proper sizes,
8461:     otherwise it will be recreated. y must be initialized to NULL if not supplied.

8463: .seealso: MatInterpolateAdd(), MatMatInterpolate(), MatMatRestrict()

8465: @*/
8466: PetscErrorCode MatMatInterpolateAdd(Mat A,Mat x,Mat w,Mat *y)
8467: {
8469:   PetscInt       M,N,Mx,Nx,Mo,My = 0,Ny = 0;
8470:   PetscBool      trans = PETSC_TRUE;
8471:   MatReuse       reuse = MAT_INITIAL_MATRIX;

8479:   MatGetSize(A,&M,&N);
8480:   MatGetSize(x,&Mx,&Nx);
8481:   if (N == Mx) trans = PETSC_FALSE;
8482:   else if (M != Mx) SETERRQ4(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Size mismatch: A %Dx%D, X %Dx%D",M,N,Mx,Nx);
8483:   Mo = trans ? N : M;
8484:   if (*y) {
8485:     MatGetSize(*y,&My,&Ny);
8486:     if (Mo == My && Nx == Ny) { reuse = MAT_REUSE_MATRIX; }
8487:     else {
8488:       if (w && *y == w) SETERRQ6(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot reuse y and w, size mismatch: A %Dx%D, X %Dx%D, Y %Dx%D",M,N,Mx,Nx,My,Ny);
8489:       MatDestroy(y);
8490:     }
8491:   }

8493:   if (w && *y == w) { /* this is to minimize changes in PCMG */
8494:     PetscBool flg;

8496:     PetscObjectQuery((PetscObject)*y,"__MatMatIntAdd_w",(PetscObject*)&w);
8497:     if (w) {
8498:       PetscInt My,Ny,Mw,Nw;

8500:       PetscObjectTypeCompare((PetscObject)*y,((PetscObject)w)->type_name,&flg);
8501:       MatGetSize(*y,&My,&Ny);
8502:       MatGetSize(w,&Mw,&Nw);
8503:       if (!flg || My != Mw || Ny != Nw) w = NULL;
8504:     }
8505:     if (!w) {
8506:       MatDuplicate(*y,MAT_COPY_VALUES,&w);
8507:       PetscObjectCompose((PetscObject)*y,"__MatMatIntAdd_w",(PetscObject)w);
8508:       PetscLogObjectParent((PetscObject)*y,(PetscObject)w);
8509:       PetscObjectDereference((PetscObject)w);
8510:     } else {
8511:       MatCopy(*y,w,UNKNOWN_NONZERO_PATTERN);
8512:     }
8513:   }
8514:   if (!trans) {
8515:     MatMatMult(A,x,reuse,PETSC_DEFAULT,y);
8516:   } else {
8517:     MatTransposeMatMult(A,x,reuse,PETSC_DEFAULT,y);
8518:   }
8519:   if (w) {
8520:     MatAXPY(*y,1.0,w,UNKNOWN_NONZERO_PATTERN);
8521:   }
8522:   return(0);
8523: }

8525: /*@
8526:    MatMatInterpolate - Y = A*X or A'*X

8528:    Neighbor-wise Collective on Mat

8530:    Input Parameters:
8531: +  mat   - the matrix
8532: -  x - the input dense matrix

8534:    Output Parameters:
8535: .  y - the output dense matrix


8538:    Level: intermediate

8540:    Notes:
8541:     This allows one to use either the restriction or interpolation (its transpose)
8542:     matrix to do the interpolation. y matrix can be reused if already created with the proper sizes,
8543:     otherwise it will be recreated. y must be initialized to NULL if not supplied.

8545: .seealso: MatInterpolate(), MatRestrict(), MatMatRestrict()

8547: @*/
8548: PetscErrorCode MatMatInterpolate(Mat A,Mat x,Mat *y)
8549: {

8553:   MatMatInterpolateAdd(A,x,NULL,y);
8554:   return(0);
8555: }

8557: /*@
8558:    MatMatRestrict - Y = A*X or A'*X

8560:    Neighbor-wise Collective on Mat

8562:    Input Parameters:
8563: +  mat   - the matrix
8564: -  x - the input dense matrix

8566:    Output Parameters:
8567: .  y - the output dense matrix


8570:    Level: intermediate

8572:    Notes:
8573:     This allows one to use either the restriction or interpolation (its transpose)
8574:     matrix to do the restriction. y matrix can be reused if already created with the proper sizes,
8575:     otherwise it will be recreated. y must be initialized to NULL if not supplied.

8577: .seealso: MatRestrict(), MatInterpolate(), MatMatInterpolate()
8578: @*/
8579: PetscErrorCode MatMatRestrict(Mat A,Mat x,Mat *y)
8580: {

8584:   MatMatInterpolateAdd(A,x,NULL,y);
8585:   return(0);
8586: }

8588: /*@
8589:    MatGetNullSpace - retrieves the null space of a matrix.

8591:    Logically Collective on Mat

8593:    Input Parameters:
8594: +  mat - the matrix
8595: -  nullsp - the null space object

8597:    Level: developer

8599: .seealso: MatCreate(), MatNullSpaceCreate(), MatSetNearNullSpace(), MatSetNullSpace()
8600: @*/
8601: PetscErrorCode MatGetNullSpace(Mat mat, MatNullSpace *nullsp)
8602: {
8606:   *nullsp = (mat->symmetric_set && mat->symmetric && !mat->nullsp) ? mat->transnullsp : mat->nullsp;
8607:   return(0);
8608: }

8610: /*@
8611:    MatSetNullSpace - attaches a null space to a matrix.

8613:    Logically Collective on Mat

8615:    Input Parameters:
8616: +  mat - the matrix
8617: -  nullsp - the null space object

8619:    Level: advanced

8621:    Notes:
8622:       This null space is used by the linear solvers. Overwrites any previous null space that may have been attached

8624:       For inconsistent singular systems (linear systems where the right hand side is not in the range of the operator) you also likely should
8625:       call MatSetTransposeNullSpace(). This allows the linear system to be solved in a least squares sense.

8627:       You can remove the null space by calling this routine with an nullsp of NULL


8630:       The fundamental theorem of linear algebra (Gilbert Strang, Introduction to Applied Mathematics, page 72) states that
8631:    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).
8632:    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
8633:    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
8634:    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).

8636:       Krylov solvers can produce the minimal norm solution to the least squares problem by utilizing MatNullSpaceRemove().

8638:     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
8639:     routine also automatically calls MatSetTransposeNullSpace().

8641: .seealso: MatCreate(), MatNullSpaceCreate(), MatSetNearNullSpace(), MatGetNullSpace(), MatSetTransposeNullSpace(), MatGetTransposeNullSpace(), MatNullSpaceRemove()
8642: @*/
8643: PetscErrorCode MatSetNullSpace(Mat mat,MatNullSpace nullsp)
8644: {

8650:   if (nullsp) {PetscObjectReference((PetscObject)nullsp);}
8651:   MatNullSpaceDestroy(&mat->nullsp);
8652:   mat->nullsp = nullsp;
8653:   if (mat->symmetric_set && mat->symmetric) {
8654:     MatSetTransposeNullSpace(mat,nullsp);
8655:   }
8656:   return(0);
8657: }

8659: /*@
8660:    MatGetTransposeNullSpace - retrieves the null space of the transpose of a matrix.

8662:    Logically Collective on Mat

8664:    Input Parameters:
8665: +  mat - the matrix
8666: -  nullsp - the null space object

8668:    Level: developer

8670: .seealso: MatCreate(), MatNullSpaceCreate(), MatSetNearNullSpace(), MatSetTransposeNullSpace(), MatSetNullSpace(), MatGetNullSpace()
8671: @*/
8672: PetscErrorCode MatGetTransposeNullSpace(Mat mat, MatNullSpace *nullsp)
8673: {
8678:   *nullsp = (mat->symmetric_set && mat->symmetric && !mat->transnullsp) ? mat->nullsp : mat->transnullsp;
8679:   return(0);
8680: }

8682: /*@
8683:    MatSetTransposeNullSpace - attaches a null space to a matrix.

8685:    Logically Collective on Mat

8687:    Input Parameters:
8688: +  mat - the matrix
8689: -  nullsp - the null space object

8691:    Level: advanced

8693:    Notes:
8694:       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.
8695:       You must also call MatSetNullSpace()


8698:       The fundamental theorem of linear algebra (Gilbert Strang, Introduction to Applied Mathematics, page 72) states that
8699:    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).
8700:    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
8701:    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
8702:    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).

8704:       Krylov solvers can produce the minimal norm solution to the least squares problem by utilizing MatNullSpaceRemove().

8706: .seealso: MatCreate(), MatNullSpaceCreate(), MatSetNearNullSpace(), MatGetNullSpace(), MatSetNullSpace(), MatGetTransposeNullSpace(), MatNullSpaceRemove()
8707: @*/
8708: PetscErrorCode MatSetTransposeNullSpace(Mat mat,MatNullSpace nullsp)
8709: {

8715:   if (nullsp) {PetscObjectReference((PetscObject)nullsp);}
8716:   MatNullSpaceDestroy(&mat->transnullsp);
8717:   mat->transnullsp = nullsp;
8718:   return(0);
8719: }

8721: /*@
8722:    MatSetNearNullSpace - attaches a null space to a matrix, which is often the null space (rigid body modes) of the operator without boundary conditions
8723:         This null space will be used to provide near null space vectors to a multigrid preconditioner built from this matrix.

8725:    Logically Collective on Mat

8727:    Input Parameters:
8728: +  mat - the matrix
8729: -  nullsp - the null space object

8731:    Level: advanced

8733:    Notes:
8734:       Overwrites any previous near null space that may have been attached

8736:       You can remove the null space by calling this routine with an nullsp of NULL

8738: .seealso: MatCreate(), MatNullSpaceCreate(), MatSetNullSpace(), MatNullSpaceCreateRigidBody(), MatGetNearNullSpace()
8739: @*/
8740: PetscErrorCode MatSetNearNullSpace(Mat mat,MatNullSpace nullsp)
8741: {

8748:   MatCheckPreallocated(mat,1);
8749:   if (nullsp) {PetscObjectReference((PetscObject)nullsp);}
8750:   MatNullSpaceDestroy(&mat->nearnullsp);
8751:   mat->nearnullsp = nullsp;
8752:   return(0);
8753: }

8755: /*@
8756:    MatGetNearNullSpace - Get null space attached with MatSetNearNullSpace()

8758:    Not Collective

8760:    Input Parameter:
8761: .  mat - the matrix

8763:    Output Parameter:
8764: .  nullsp - the null space object, NULL if not set

8766:    Level: developer

8768: .seealso: MatSetNearNullSpace(), MatGetNullSpace(), MatNullSpaceCreate()
8769: @*/
8770: PetscErrorCode MatGetNearNullSpace(Mat mat,MatNullSpace *nullsp)
8771: {
8776:   MatCheckPreallocated(mat,1);
8777:   *nullsp = mat->nearnullsp;
8778:   return(0);
8779: }

8781: /*@C
8782:    MatICCFactor - Performs in-place incomplete Cholesky factorization of matrix.

8784:    Collective on Mat

8786:    Input Parameters:
8787: +  mat - the matrix
8788: .  row - row/column permutation
8789: .  fill - expected fill factor >= 1.0
8790: -  level - level of fill, for ICC(k)

8792:    Notes:
8793:    Probably really in-place only when level of fill is zero, otherwise allocates
8794:    new space to store factored matrix and deletes previous memory.

8796:    Most users should employ the simplified KSP interface for linear solvers
8797:    instead of working directly with matrix algebra routines such as this.
8798:    See, e.g., KSPCreate().

8800:    Level: developer


8803: .seealso: MatICCFactorSymbolic(), MatLUFactorNumeric(), MatCholeskyFactor()

8805:     Developer Note: fortran interface is not autogenerated as the f90
8806:     interface defintion cannot be generated correctly [due to MatFactorInfo]

8808: @*/
8809: PetscErrorCode MatICCFactor(Mat mat,IS row,const MatFactorInfo *info)
8810: {

8818:   if (mat->rmap->N != mat->cmap->N) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONG,"matrix must be square");
8819:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
8820:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
8821:   if (!mat->ops->iccfactor) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
8822:   MatCheckPreallocated(mat,1);
8823:   (*mat->ops->iccfactor)(mat,row,info);
8824:   PetscObjectStateIncrease((PetscObject)mat);
8825:   return(0);
8826: }

8828: /*@
8829:    MatDiagonalScaleLocal - Scales columns of a matrix given the scaling values including the
8830:          ghosted ones.

8832:    Not Collective

8834:    Input Parameters:
8835: +  mat - the matrix
8836: -  diag = the diagonal values, including ghost ones

8838:    Level: developer

8840:    Notes:
8841:     Works only for MPIAIJ and MPIBAIJ matrices

8843: .seealso: MatDiagonalScale()
8844: @*/
8845: PetscErrorCode MatDiagonalScaleLocal(Mat mat,Vec diag)
8846: {
8848:   PetscMPIInt    size;


8855:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Matrix must be already assembled");
8856:   PetscLogEventBegin(MAT_Scale,mat,0,0,0);
8857:   MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);
8858:   if (size == 1) {
8859:     PetscInt n,m;
8860:     VecGetSize(diag,&n);
8861:     MatGetSize(mat,NULL,&m);
8862:     if (m == n) {
8863:       MatDiagonalScale(mat,NULL,diag);
8864:     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only supported for sequential matrices when no ghost points/periodic conditions");
8865:   } else {
8866:     PetscUseMethod(mat,"MatDiagonalScaleLocal_C",(Mat,Vec),(mat,diag));
8867:   }
8868:   PetscLogEventEnd(MAT_Scale,mat,0,0,0);
8869:   PetscObjectStateIncrease((PetscObject)mat);
8870:   return(0);
8871: }

8873: /*@
8874:    MatGetInertia - Gets the inertia from a factored matrix

8876:    Collective on Mat

8878:    Input Parameter:
8879: .  mat - the matrix

8881:    Output Parameters:
8882: +   nneg - number of negative eigenvalues
8883: .   nzero - number of zero eigenvalues
8884: -   npos - number of positive eigenvalues

8886:    Level: advanced

8888:    Notes:
8889:     Matrix must have been factored by MatCholeskyFactor()


8892: @*/
8893: PetscErrorCode MatGetInertia(Mat mat,PetscInt *nneg,PetscInt *nzero,PetscInt *npos)
8894: {

8900:   if (!mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Unfactored matrix");
8901:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Numeric factor mat is not assembled");
8902:   if (!mat->ops->getinertia) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
8903:   (*mat->ops->getinertia)(mat,nneg,nzero,npos);
8904:   return(0);
8905: }

8907: /* ----------------------------------------------------------------*/
8908: /*@C
8909:    MatSolves - Solves A x = b, given a factored matrix, for a collection of vectors

8911:    Neighbor-wise Collective on Mats

8913:    Input Parameters:
8914: +  mat - the factored matrix
8915: -  b - the right-hand-side vectors

8917:    Output Parameter:
8918: .  x - the result vectors

8920:    Notes:
8921:    The vectors b and x cannot be the same.  I.e., one cannot
8922:    call MatSolves(A,x,x).

8924:    Notes:
8925:    Most users should employ the simplified KSP interface for linear solvers
8926:    instead of working directly with matrix algebra routines such as this.
8927:    See, e.g., KSPCreate().

8929:    Level: developer

8931: .seealso: MatSolveAdd(), MatSolveTranspose(), MatSolveTransposeAdd(), MatSolve()
8932: @*/
8933: PetscErrorCode MatSolves(Mat mat,Vecs b,Vecs x)
8934: {

8940:   if (x == b) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_IDN,"x and b must be different vectors");
8941:   if (!mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Unfactored matrix");
8942:   if (!mat->rmap->N && !mat->cmap->N) return(0);

8944:   if (!mat->ops->solves) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)mat)->type_name);
8945:   MatCheckPreallocated(mat,1);
8946:   PetscLogEventBegin(MAT_Solves,mat,0,0,0);
8947:   (*mat->ops->solves)(mat,b,x);
8948:   PetscLogEventEnd(MAT_Solves,mat,0,0,0);
8949:   return(0);
8950: }

8952: /*@
8953:    MatIsSymmetric - Test whether a matrix is symmetric

8955:    Collective on Mat

8957:    Input Parameter:
8958: +  A - the matrix to test
8959: -  tol - difference between value and its transpose less than this amount counts as equal (use 0.0 for exact transpose)

8961:    Output Parameters:
8962: .  flg - the result

8964:    Notes:
8965:     For real numbers MatIsSymmetric() and MatIsHermitian() return identical results

8967:    Level: intermediate

8969: .seealso: MatTranspose(), MatIsTranspose(), MatIsHermitian(), MatIsStructurallySymmetric(), MatSetOption(), MatIsSymmetricKnown()
8970: @*/
8971: PetscErrorCode MatIsSymmetric(Mat A,PetscReal tol,PetscBool  *flg)
8972: {


8979:   if (!A->symmetric_set) {
8980:     if (!A->ops->issymmetric) {
8981:       MatType mattype;
8982:       MatGetType(A,&mattype);
8983:       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix of type %s does not support checking for symmetric",mattype);
8984:     }
8985:     (*A->ops->issymmetric)(A,tol,flg);
8986:     if (!tol) {
8987:       MatSetOption(A,MAT_SYMMETRIC,*flg);
8988:     }
8989:   } else if (A->symmetric) {
8990:     *flg = PETSC_TRUE;
8991:   } else if (!tol) {
8992:     *flg = PETSC_FALSE;
8993:   } else {
8994:     if (!A->ops->issymmetric) {
8995:       MatType mattype;
8996:       MatGetType(A,&mattype);
8997:       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix of type %s does not support checking for symmetric",mattype);
8998:     }
8999:     (*A->ops->issymmetric)(A,tol,flg);
9000:   }
9001:   return(0);
9002: }

9004: /*@
9005:    MatIsHermitian - Test whether a matrix is Hermitian

9007:    Collective on Mat

9009:    Input Parameter:
9010: +  A - the matrix to test
9011: -  tol - difference between value and its transpose less than this amount counts as equal (use 0.0 for exact Hermitian)

9013:    Output Parameters:
9014: .  flg - the result

9016:    Level: intermediate

9018: .seealso: MatTranspose(), MatIsTranspose(), MatIsHermitian(), MatIsStructurallySymmetric(), MatSetOption(),
9019:           MatIsSymmetricKnown(), MatIsSymmetric()
9020: @*/
9021: PetscErrorCode MatIsHermitian(Mat A,PetscReal tol,PetscBool  *flg)
9022: {


9029:   if (!A->hermitian_set) {
9030:     if (!A->ops->ishermitian) {
9031:       MatType mattype;
9032:       MatGetType(A,&mattype);
9033:       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix of type %s does not support checking for hermitian",mattype);
9034:     }
9035:     (*A->ops->ishermitian)(A,tol,flg);
9036:     if (!tol) {
9037:       MatSetOption(A,MAT_HERMITIAN,*flg);
9038:     }
9039:   } else if (A->hermitian) {
9040:     *flg = PETSC_TRUE;
9041:   } else if (!tol) {
9042:     *flg = PETSC_FALSE;
9043:   } else {
9044:     if (!A->ops->ishermitian) {
9045:       MatType mattype;
9046:       MatGetType(A,&mattype);
9047:       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix of type %s does not support checking for hermitian",mattype);
9048:     }
9049:     (*A->ops->ishermitian)(A,tol,flg);
9050:   }
9051:   return(0);
9052: }

9054: /*@
9055:    MatIsSymmetricKnown - Checks the flag on the matrix to see if it is symmetric.

9057:    Not Collective

9059:    Input Parameter:
9060: .  A - the matrix to check

9062:    Output Parameters:
9063: +  set - if the symmetric flag is set (this tells you if the next flag is valid)
9064: -  flg - the result

9066:    Level: advanced

9068:    Note: Does not check the matrix values directly, so this may return unknown (set = PETSC_FALSE). Use MatIsSymmetric()
9069:          if you want it explicitly checked

9071: .seealso: MatTranspose(), MatIsTranspose(), MatIsHermitian(), MatIsStructurallySymmetric(), MatSetOption(), MatIsSymmetric()
9072: @*/
9073: PetscErrorCode MatIsSymmetricKnown(Mat A,PetscBool *set,PetscBool *flg)
9074: {
9079:   if (A->symmetric_set) {
9080:     *set = PETSC_TRUE;
9081:     *flg = A->symmetric;
9082:   } else {
9083:     *set = PETSC_FALSE;
9084:   }
9085:   return(0);
9086: }

9088: /*@
9089:    MatIsHermitianKnown - Checks the flag on the matrix to see if it is hermitian.

9091:    Not Collective

9093:    Input Parameter:
9094: .  A - the matrix to check

9096:    Output Parameters:
9097: +  set - if the hermitian flag is set (this tells you if the next flag is valid)
9098: -  flg - the result

9100:    Level: advanced

9102:    Note: Does not check the matrix values directly, so this may return unknown (set = PETSC_FALSE). Use MatIsHermitian()
9103:          if you want it explicitly checked

9105: .seealso: MatTranspose(), MatIsTranspose(), MatIsHermitian(), MatIsStructurallySymmetric(), MatSetOption(), MatIsSymmetric()
9106: @*/
9107: PetscErrorCode MatIsHermitianKnown(Mat A,PetscBool *set,PetscBool *flg)
9108: {
9113:   if (A->hermitian_set) {
9114:     *set = PETSC_TRUE;
9115:     *flg = A->hermitian;
9116:   } else {
9117:     *set = PETSC_FALSE;
9118:   }
9119:   return(0);
9120: }

9122: /*@
9123:    MatIsStructurallySymmetric - Test whether a matrix is structurally symmetric

9125:    Collective on Mat

9127:    Input Parameter:
9128: .  A - the matrix to test

9130:    Output Parameters:
9131: .  flg - the result

9133:    Level: intermediate

9135: .seealso: MatTranspose(), MatIsTranspose(), MatIsHermitian(), MatIsSymmetric(), MatSetOption()
9136: @*/
9137: PetscErrorCode MatIsStructurallySymmetric(Mat A,PetscBool *flg)
9138: {

9144:   if (!A->structurally_symmetric_set) {
9145:     if (!A->ops->isstructurallysymmetric) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Matrix of type %s does not support checking for structural symmetric",((PetscObject)A)->type_name);
9146:     (*A->ops->isstructurallysymmetric)(A,flg);
9147:     MatSetOption(A,MAT_STRUCTURALLY_SYMMETRIC,*flg);
9148:   } else *flg = A->structurally_symmetric;
9149:   return(0);
9150: }

9152: /*@
9153:    MatStashGetInfo - Gets how many values are currently in the matrix stash, i.e. need
9154:        to be communicated to other processors during the MatAssemblyBegin/End() process

9156:     Not collective

9158:    Input Parameter:
9159: .   vec - the vector

9161:    Output Parameters:
9162: +   nstash   - the size of the stash
9163: .   reallocs - the number of additional mallocs incurred.
9164: .   bnstash   - the size of the block stash
9165: -   breallocs - the number of additional mallocs incurred.in the block stash

9167:    Level: advanced

9169: .seealso: MatAssemblyBegin(), MatAssemblyEnd(), Mat, MatStashSetInitialSize()

9171: @*/
9172: PetscErrorCode MatStashGetInfo(Mat mat,PetscInt *nstash,PetscInt *reallocs,PetscInt *bnstash,PetscInt *breallocs)
9173: {

9177:   MatStashGetInfo_Private(&mat->stash,nstash,reallocs);
9178:   MatStashGetInfo_Private(&mat->bstash,bnstash,breallocs);
9179:   return(0);
9180: }

9182: /*@C
9183:    MatCreateVecs - Get vector(s) compatible with the matrix, i.e. with the same
9184:      parallel layout

9186:    Collective on Mat

9188:    Input Parameter:
9189: .  mat - the matrix

9191:    Output Parameter:
9192: +   right - (optional) vector that the matrix can be multiplied against
9193: -   left - (optional) vector that the matrix vector product can be stored in

9195:    Notes:
9196:     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().

9198:   Notes:
9199:     These are new vectors which are not owned by the Mat, they should be destroyed in VecDestroy() when no longer needed

9201:   Level: advanced

9203: .seealso: MatCreate(), VecDestroy()
9204: @*/
9205: PetscErrorCode MatCreateVecs(Mat mat,Vec *right,Vec *left)
9206: {

9212:   if (mat->ops->getvecs) {
9213:     (*mat->ops->getvecs)(mat,right,left);
9214:   } else {
9215:     PetscInt rbs,cbs;
9216:     MatGetBlockSizes(mat,&rbs,&cbs);
9217:     if (right) {
9218:       if (mat->cmap->n < 0) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"PetscLayout for columns not yet setup");
9219:       VecCreate(PetscObjectComm((PetscObject)mat),right);
9220:       VecSetSizes(*right,mat->cmap->n,PETSC_DETERMINE);
9221:       VecSetBlockSize(*right,cbs);
9222:       VecSetType(*right,mat->defaultvectype);
9223:       PetscLayoutReference(mat->cmap,&(*right)->map);
9224:     }
9225:     if (left) {
9226:       if (mat->rmap->n < 0) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"PetscLayout for rows not yet setup");
9227:       VecCreate(PetscObjectComm((PetscObject)mat),left);
9228:       VecSetSizes(*left,mat->rmap->n,PETSC_DETERMINE);
9229:       VecSetBlockSize(*left,rbs);
9230:       VecSetType(*left,mat->defaultvectype);
9231:       PetscLayoutReference(mat->rmap,&(*left)->map);
9232:     }
9233:   }
9234:   return(0);
9235: }

9237: /*@C
9238:    MatFactorInfoInitialize - Initializes a MatFactorInfo data structure
9239:      with default values.

9241:    Not Collective

9243:    Input Parameters:
9244: .    info - the MatFactorInfo data structure


9247:    Notes:
9248:     The solvers are generally used through the KSP and PC objects, for example
9249:           PCLU, PCILU, PCCHOLESKY, PCICC

9251:    Level: developer

9253: .seealso: MatFactorInfo

9255:     Developer Note: fortran interface is not autogenerated as the f90
9256:     interface defintion cannot be generated correctly [due to MatFactorInfo]

9258: @*/

9260: PetscErrorCode MatFactorInfoInitialize(MatFactorInfo *info)
9261: {

9265:   PetscMemzero(info,sizeof(MatFactorInfo));
9266:   return(0);
9267: }

9269: /*@
9270:    MatFactorSetSchurIS - Set indices corresponding to the Schur complement you wish to have computed

9272:    Collective on Mat

9274:    Input Parameters:
9275: +  mat - the factored matrix
9276: -  is - the index set defining the Schur indices (0-based)

9278:    Notes:
9279:     Call MatFactorSolveSchurComplement() or MatFactorSolveSchurComplementTranspose() after this call to solve a Schur complement system.

9281:    You can call MatFactorGetSchurComplement() or MatFactorCreateSchurComplement() after this call.

9283:    Level: developer

9285: .seealso: MatGetFactor(), MatFactorGetSchurComplement(), MatFactorRestoreSchurComplement(), MatFactorCreateSchurComplement(), MatFactorSolveSchurComplement(),
9286:           MatFactorSolveSchurComplementTranspose(), MatFactorSolveSchurComplement()

9288: @*/
9289: PetscErrorCode MatFactorSetSchurIS(Mat mat,IS is)
9290: {
9291:   PetscErrorCode ierr,(*f)(Mat,IS);

9299:   if (!mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
9300:   PetscObjectQueryFunction((PetscObject)mat,"MatFactorSetSchurIS_C",&f);
9301:   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");
9302:   MatDestroy(&mat->schur);
9303:   (*f)(mat,is);
9304:   if (!mat->schur) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_PLIB,"Schur complement has not been created");
9305:   return(0);
9306: }

9308: /*@
9309:   MatFactorCreateSchurComplement - Create a Schur complement matrix object using Schur data computed during the factorization step

9311:    Logically Collective on Mat

9313:    Input Parameters:
9314: +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
9315: .  S - location where to return the Schur complement, can be NULL
9316: -  status - the status of the Schur complement matrix, can be NULL

9318:    Notes:
9319:    You must call MatFactorSetSchurIS() before calling this routine.

9321:    The routine provides a copy of the Schur matrix stored within the solver data structures.
9322:    The caller must destroy the object when it is no longer needed.
9323:    If MatFactorInvertSchurComplement() has been called, the routine gets back the inverse.

9325:    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)

9327:    Developer Notes:
9328:     The reason this routine exists is because the representation of the Schur complement within the factor matrix may be different than a standard PETSc
9329:    matrix representation and we normally do not want to use the time or memory to make a copy as a regular PETSc matrix.

9331:    See MatCreateSchurComplement() or MatGetSchurComplement() for ways to create virtual or approximate Schur complements.

9333:    Level: advanced

9335:    References:

9337: .seealso: MatGetFactor(), MatFactorSetSchurIS(), MatFactorGetSchurComplement(), MatFactorSchurStatus
9338: @*/
9339: PetscErrorCode MatFactorCreateSchurComplement(Mat F,Mat* S,MatFactorSchurStatus* status)
9340: {

9347:   if (S) {
9348:     PetscErrorCode (*f)(Mat,Mat*);

9350:     PetscObjectQueryFunction((PetscObject)F,"MatFactorCreateSchurComplement_C",&f);
9351:     if (f) {
9352:       (*f)(F,S);
9353:     } else {
9354:       MatDuplicate(F->schur,MAT_COPY_VALUES,S);
9355:     }
9356:   }
9357:   if (status) *status = F->schur_status;
9358:   return(0);
9359: }

9361: /*@
9362:   MatFactorGetSchurComplement - Gets access to a Schur complement matrix using the current Schur data within a factored matrix

9364:    Logically Collective on Mat

9366:    Input Parameters:
9367: +  F - the factored matrix obtained by calling MatGetFactor()
9368: .  *S - location where to return the Schur complement, can be NULL
9369: -  status - the status of the Schur complement matrix, can be NULL

9371:    Notes:
9372:    You must call MatFactorSetSchurIS() before calling this routine.

9374:    Schur complement mode is currently implemented for sequential matrices.
9375:    The routine returns a the Schur Complement stored within the data strutures of the solver.
9376:    If MatFactorInvertSchurComplement() has previously been called, the returned matrix is actually the inverse of the Schur complement.
9377:    The returned matrix should not be destroyed; the caller should call MatFactorRestoreSchurComplement() when the object is no longer needed.

9379:    Use MatFactorCreateSchurComplement() to create a copy of the Schur complement matrix that is within a factored matrix

9381:    See MatCreateSchurComplement() or MatGetSchurComplement() for ways to create virtual or approximate Schur complements.

9383:    Level: advanced

9385:    References:

9387: .seealso: MatGetFactor(), MatFactorSetSchurIS(), MatFactorRestoreSchurComplement(), MatFactorCreateSchurComplement(), MatFactorSchurStatus
9388: @*/
9389: PetscErrorCode MatFactorGetSchurComplement(Mat F,Mat* S,MatFactorSchurStatus* status)
9390: {
9395:   if (S) *S = F->schur;
9396:   if (status) *status = F->schur_status;
9397:   return(0);
9398: }

9400: /*@
9401:   MatFactorRestoreSchurComplement - Restore the Schur complement matrix object obtained from a call to MatFactorGetSchurComplement

9403:    Logically Collective on Mat

9405:    Input Parameters:
9406: +  F - the factored matrix obtained by calling MatGetFactor()
9407: .  *S - location where the Schur complement is stored
9408: -  status - the status of the Schur complement matrix (see MatFactorSchurStatus)

9410:    Notes:

9412:    Level: advanced

9414:    References:

9416: .seealso: MatGetFactor(), MatFactorSetSchurIS(), MatFactorRestoreSchurComplement(), MatFactorCreateSchurComplement(), MatFactorSchurStatus
9417: @*/
9418: PetscErrorCode MatFactorRestoreSchurComplement(Mat F,Mat* S,MatFactorSchurStatus status)
9419: {

9424:   if (S) {
9426:     *S = NULL;
9427:   }
9428:   F->schur_status = status;
9429:   MatFactorUpdateSchurStatus_Private(F);
9430:   return(0);
9431: }

9433: /*@
9434:   MatFactorSolveSchurComplementTranspose - Solve the transpose of the Schur complement system computed during the factorization step

9436:    Logically Collective on Mat

9438:    Input Parameters:
9439: +  F - the factored matrix obtained by calling MatGetFactor()
9440: .  rhs - location where the right hand side of the Schur complement system is stored
9441: -  sol - location where the solution of the Schur complement system has to be returned

9443:    Notes:
9444:    The sizes of the vectors should match the size of the Schur complement

9446:    Must be called after MatFactorSetSchurIS()

9448:    Level: advanced

9450:    References:

9452: .seealso: MatGetFactor(), MatFactorSetSchurIS(), MatFactorSolveSchurComplement()
9453: @*/
9454: PetscErrorCode MatFactorSolveSchurComplementTranspose(Mat F, Vec rhs, Vec sol)
9455: {

9467:   MatFactorFactorizeSchurComplement(F);
9468:   switch (F->schur_status) {
9469:   case MAT_FACTOR_SCHUR_FACTORED:
9470:     MatSolveTranspose(F->schur,rhs,sol);
9471:     break;
9472:   case MAT_FACTOR_SCHUR_INVERTED:
9473:     MatMultTranspose(F->schur,rhs,sol);
9474:     break;
9475:   default:
9476:     SETERRQ1(PetscObjectComm((PetscObject)F),PETSC_ERR_SUP,"Unhandled MatFactorSchurStatus %D",F->schur_status);
9477:   }
9478:   return(0);
9479: }

9481: /*@
9482:   MatFactorSolveSchurComplement - Solve the Schur complement system computed during the factorization step

9484:    Logically Collective on Mat

9486:    Input Parameters:
9487: +  F - the factored matrix obtained by calling MatGetFactor()
9488: .  rhs - location where the right hand side of the Schur complement system is stored
9489: -  sol - location where the solution of the Schur complement system has to be returned

9491:    Notes:
9492:    The sizes of the vectors should match the size of the Schur complement

9494:    Must be called after MatFactorSetSchurIS()

9496:    Level: advanced

9498:    References:

9500: .seealso: MatGetFactor(), MatFactorSetSchurIS(), MatFactorSolveSchurComplementTranspose()
9501: @*/
9502: PetscErrorCode MatFactorSolveSchurComplement(Mat F, Vec rhs, Vec sol)
9503: {

9515:   MatFactorFactorizeSchurComplement(F);
9516:   switch (F->schur_status) {
9517:   case MAT_FACTOR_SCHUR_FACTORED:
9518:     MatSolve(F->schur,rhs,sol);
9519:     break;
9520:   case MAT_FACTOR_SCHUR_INVERTED:
9521:     MatMult(F->schur,rhs,sol);
9522:     break;
9523:   default:
9524:     SETERRQ1(PetscObjectComm((PetscObject)F),PETSC_ERR_SUP,"Unhandled MatFactorSchurStatus %D",F->schur_status);
9525:   }
9526:   return(0);
9527: }

9529: /*@
9530:   MatFactorInvertSchurComplement - Invert the Schur complement matrix computed during the factorization step

9532:    Logically Collective on Mat

9534:    Input Parameters:
9535: .  F - the factored matrix obtained by calling MatGetFactor()

9537:    Notes:
9538:     Must be called after MatFactorSetSchurIS().

9540:    Call MatFactorGetSchurComplement() or  MatFactorCreateSchurComplement() AFTER this call to actually compute the inverse and get access to it.

9542:    Level: advanced

9544:    References:

9546: .seealso: MatGetFactor(), MatFactorSetSchurIS(), MatFactorGetSchurComplement(), MatFactorCreateSchurComplement()
9547: @*/
9548: PetscErrorCode MatFactorInvertSchurComplement(Mat F)
9549: {

9555:   if (F->schur_status == MAT_FACTOR_SCHUR_INVERTED) return(0);
9556:   MatFactorFactorizeSchurComplement(F);
9557:   MatFactorInvertSchurComplement_Private(F);
9558:   F->schur_status = MAT_FACTOR_SCHUR_INVERTED;
9559:   return(0);
9560: }

9562: /*@
9563:   MatFactorFactorizeSchurComplement - Factorize the Schur complement matrix computed during the factorization step

9565:    Logically Collective on Mat

9567:    Input Parameters:
9568: .  F - the factored matrix obtained by calling MatGetFactor()

9570:    Notes:
9571:     Must be called after MatFactorSetSchurIS().

9573:    Level: advanced

9575:    References:

9577: .seealso: MatGetFactor(), MatFactorSetSchurIS(), MatFactorInvertSchurComplement()
9578: @*/
9579: PetscErrorCode MatFactorFactorizeSchurComplement(Mat F)
9580: {

9586:   if (F->schur_status == MAT_FACTOR_SCHUR_INVERTED || F->schur_status == MAT_FACTOR_SCHUR_FACTORED) return(0);
9587:   MatFactorFactorizeSchurComplement_Private(F);
9588:   F->schur_status = MAT_FACTOR_SCHUR_FACTORED;
9589:   return(0);
9590: }

9592: /*@
9593:    MatPtAP - Creates the matrix product C = P^T * A * P

9595:    Neighbor-wise Collective on Mat

9597:    Input Parameters:
9598: +  A - the matrix
9599: .  P - the projection matrix
9600: .  scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
9601: -  fill - expected fill as ratio of nnz(C)/(nnz(A) + nnz(P)), use PETSC_DEFAULT if you do not have a good estimate
9602:           if the result is a dense matrix this is irrelevent

9604:    Output Parameters:
9605: .  C - the product matrix

9607:    Notes:
9608:    C will be created and must be destroyed by the user with MatDestroy().

9610:    For matrix types without special implementation the function fallbacks to MatMatMult() followed by MatTransposeMatMult().

9612:    Level: intermediate

9614: .seealso: MatMatMult(), MatRARt()
9615: @*/
9616: PetscErrorCode MatPtAP(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
9617: {

9621:   if (scall == MAT_REUSE_MATRIX) MatCheckProduct(*C,5);
9622:   if (scall == MAT_INPLACE_MATRIX) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Inplace product not supported");

9624:   if (scall == MAT_INITIAL_MATRIX) {
9625:     MatProductCreate(A,P,NULL,C);
9626:     MatProductSetType(*C,MATPRODUCT_PtAP);
9627:     MatProductSetAlgorithm(*C,"default");
9628:     MatProductSetFill(*C,fill);

9630:     (*C)->product->api_user = PETSC_TRUE;
9631:     MatProductSetFromOptions(*C);
9632:     if (!(*C)->ops->productsymbolic) SETERRQ3(PetscObjectComm((PetscObject)(*C)),PETSC_ERR_SUP,"MatProduct %s not supported for A %s and P %s",MatProductTypes[MATPRODUCT_PtAP],((PetscObject)A)->type_name,((PetscObject)P)->type_name);
9633:     MatProductSymbolic(*C);
9634:   } else { /* scall == MAT_REUSE_MATRIX */
9635:     MatProductReplaceMats(A,P,NULL,*C);
9636:   }

9638:   MatProductNumeric(*C);
9639:   if (A->symmetric_set && A->symmetric) {
9640:     MatSetOption(*C,MAT_SYMMETRIC,PETSC_TRUE);
9641:   }
9642:   return(0);
9643: }

9645: /*@
9646:    MatRARt - Creates the matrix product C = R * A * R^T

9648:    Neighbor-wise Collective on Mat

9650:    Input Parameters:
9651: +  A - the matrix
9652: .  R - the projection matrix
9653: .  scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
9654: -  fill - expected fill as ratio of nnz(C)/nnz(A), use PETSC_DEFAULT if you do not have a good estimate
9655:           if the result is a dense matrix this is irrelevent

9657:    Output Parameters:
9658: .  C - the product matrix

9660:    Notes:
9661:    C will be created and must be destroyed by the user with MatDestroy().

9663:    This routine is currently only implemented for pairs of AIJ matrices and classes
9664:    which inherit from AIJ. Due to PETSc sparse matrix block row distribution among processes,
9665:    parallel MatRARt is implemented via explicit transpose of R, which could be very expensive.
9666:    We recommend using MatPtAP().

9668:    Level: intermediate

9670: .seealso: MatMatMult(), MatPtAP()
9671: @*/
9672: PetscErrorCode MatRARt(Mat A,Mat R,MatReuse scall,PetscReal fill,Mat *C)
9673: {

9677:   if (scall == MAT_REUSE_MATRIX) MatCheckProduct(*C,5);
9678:   if (scall == MAT_INPLACE_MATRIX) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Inplace product not supported");

9680:   if (scall == MAT_INITIAL_MATRIX) {
9681:     MatProductCreate(A,R,NULL,C);
9682:     MatProductSetType(*C,MATPRODUCT_RARt);
9683:     MatProductSetAlgorithm(*C,"default");
9684:     MatProductSetFill(*C,fill);

9686:     (*C)->product->api_user = PETSC_TRUE;
9687:     MatProductSetFromOptions(*C);
9688:     if (!(*C)->ops->productsymbolic) SETERRQ3(PetscObjectComm((PetscObject)(*C)),PETSC_ERR_SUP,"MatProduct %s not supported for A %s and R %s",MatProductTypes[MATPRODUCT_RARt],((PetscObject)A)->type_name,((PetscObject)R)->type_name);
9689:     MatProductSymbolic(*C);
9690:   } else { /* scall == MAT_REUSE_MATRIX */
9691:     MatProductReplaceMats(A,R,NULL,*C);
9692:   }

9694:   MatProductNumeric(*C);
9695:   if (A->symmetric_set && A->symmetric) {
9696:     MatSetOption(*C,MAT_SYMMETRIC,PETSC_TRUE);
9697:   }
9698:   return(0);
9699: }


9702: static PetscErrorCode MatProduct_Private(Mat A,Mat B,MatReuse scall,PetscReal fill,MatProductType ptype, Mat *C)
9703: {

9707:   if (scall == MAT_INPLACE_MATRIX) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Inplace product not supported");

9709:   if (scall == MAT_INITIAL_MATRIX) {
9710:     PetscInfo1(A,"Calling MatProduct API with MAT_INITIAL_MATRIX and product type %s\n",MatProductTypes[ptype]);
9711:     MatProductCreate(A,B,NULL,C);
9712:     MatProductSetType(*C,ptype);
9713:     MatProductSetAlgorithm(*C,MATPRODUCTALGORITHM_DEFAULT);
9714:     MatProductSetFill(*C,fill);

9716:     (*C)->product->api_user = PETSC_TRUE;
9717:     MatProductSetFromOptions(*C);
9718:     MatProductSymbolic(*C);
9719:   } else { /* scall == MAT_REUSE_MATRIX */
9720:     Mat_Product *product = (*C)->product;

9722:     PetscInfo2(A,"Calling MatProduct API with MAT_REUSE_MATRIX %s product present and product type %s\n",product ? "with" : "without",MatProductTypes[ptype]);
9723:     if (!product) {
9724:       /* user provide the dense matrix *C without calling MatProductCreate() */
9725:       PetscBool isdense;

9727:       PetscObjectBaseTypeCompareAny((PetscObject)(*C),&isdense,MATSEQDENSE,MATMPIDENSE,"");
9728:       if (isdense) {
9729:         /* user wants to reuse an assembled dense matrix */
9730:         /* Create product -- see MatCreateProduct() */
9731:         MatProductCreate_Private(A,B,NULL,*C);
9732:         product = (*C)->product;
9733:         product->fill     = fill;
9734:         product->api_user = PETSC_TRUE;
9735:         product->clear    = PETSC_TRUE;

9737:         MatProductSetType(*C,ptype);
9738:         MatProductSetFromOptions(*C);
9739:         if (!(*C)->ops->productsymbolic) SETERRQ3(PetscObjectComm((PetscObject)(*C)),PETSC_ERR_SUP,"MatProduct %s not supported for %s and %s",MatProductTypes[ptype],((PetscObject)A)->type_name,((PetscObject)B)->type_name);
9740:         MatProductSymbolic(*C);
9741:       } else SETERRQ(PetscObjectComm((PetscObject)(*C)),PETSC_ERR_SUP,"Call MatProductCreate() first");
9742:     } else { /* user may change input matrices A or B when REUSE */
9743:       MatProductReplaceMats(A,B,NULL,*C);
9744:     }
9745:   }
9746:   MatProductNumeric(*C);
9747:   return(0);
9748: }

9750: /*@
9751:    MatMatMult - Performs Matrix-Matrix Multiplication C=A*B.

9753:    Neighbor-wise Collective on Mat

9755:    Input Parameters:
9756: +  A - the left matrix
9757: .  B - the right matrix
9758: .  scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
9759: -  fill - expected fill as ratio of nnz(C)/(nnz(A) + nnz(B)), use PETSC_DEFAULT if you do not have a good estimate
9760:           if the result is a dense matrix this is irrelevent

9762:    Output Parameters:
9763: .  C - the product matrix

9765:    Notes:
9766:    Unless scall is MAT_REUSE_MATRIX C will be created.

9768:    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
9769:    call to this function with MAT_INITIAL_MATRIX.

9771:    To determine the correct fill value, run with -info and search for the string "Fill ratio" to see the value actually needed.

9773:    If you have many matrices with the same non-zero structure to multiply, you should use MatProductCreate()/MatProductSymbolic(C)/ReplaceMats(), and call MatProductNumeric() repeatedly.

9775:    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 with MAT_REUSE_MATRIX, rather than first having MatMatMult() create it for you. You can NEVER do this if the matrix C is sparse.

9777:    Level: intermediate

9779: .seealso: MatTransposeMatMult(), MatMatTransposeMult(), MatPtAP()
9780: @*/
9781: PetscErrorCode MatMatMult(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
9782: {

9786:   MatProduct_Private(A,B,scall,fill,MATPRODUCT_AB,C);
9787:   return(0);
9788: }

9790: /*@
9791:    MatMatTransposeMult - Performs Matrix-Matrix Multiplication C=A*B^T.

9793:    Neighbor-wise Collective on Mat

9795:    Input Parameters:
9796: +  A - the left matrix
9797: .  B - the right matrix
9798: .  scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
9799: -  fill - expected fill as ratio of nnz(C)/(nnz(A) + nnz(B)), use PETSC_DEFAULT if not known

9801:    Output Parameters:
9802: .  C - the product matrix

9804:    Notes:
9805:    C will be created if MAT_INITIAL_MATRIX and must be destroyed by the user with MatDestroy().

9807:    MAT_REUSE_MATRIX can only be used if the matrices A and B have the same nonzero pattern as in the previous call

9809:   To determine the correct fill value, run with -info and search for the string "Fill ratio" to see the value
9810:    actually needed.

9812:    This routine is currently only implemented for pairs of SeqAIJ matrices, for the SeqDense class,
9813:    and for pairs of MPIDense matrices.

9815:    Options Database Keys:
9816: .  -matmattransmult_mpidense_mpidense_via {allgatherv,cyclic} - Choose between algorthims for MPIDense matrices: the
9817:                                                                 first redundantly copies the transposed B matrix on each process and requiers O(log P) communication complexity;
9818:                                                                 the second never stores more than one portion of the B matrix at a time by requires O(P) communication complexity.

9820:    Level: intermediate

9822: .seealso: MatMatMult(), MatTransposeMatMult() MatPtAP()
9823: @*/
9824: PetscErrorCode MatMatTransposeMult(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
9825: {

9829:   MatProduct_Private(A,B,scall,fill,MATPRODUCT_ABt,C);
9830:   return(0);
9831: }

9833: /*@
9834:    MatTransposeMatMult - Performs Matrix-Matrix Multiplication C=A^T*B.

9836:    Neighbor-wise Collective on Mat

9838:    Input Parameters:
9839: +  A - the left matrix
9840: .  B - the right matrix
9841: .  scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
9842: -  fill - expected fill as ratio of nnz(C)/(nnz(A) + nnz(B)), use PETSC_DEFAULT if not known

9844:    Output Parameters:
9845: .  C - the product matrix

9847:    Notes:
9848:    C will be created if MAT_INITIAL_MATRIX and must be destroyed by the user with MatDestroy().

9850:    MAT_REUSE_MATRIX can only be used if the matrices A and B have the same nonzero pattern as in the previous call.

9852:   To determine the correct fill value, run with -info and search for the string "Fill ratio" to see the value
9853:    actually needed.

9855:    This routine is currently implemented for pairs of AIJ matrices and pairs of SeqDense matrices and classes
9856:    which inherit from SeqAIJ.  C will be of same type as the input matrices.

9858:    Level: intermediate

9860: .seealso: MatMatMult(), MatMatTransposeMult(), MatPtAP()
9861: @*/
9862: PetscErrorCode MatTransposeMatMult(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
9863: {

9867:   MatProduct_Private(A,B,scall,fill,MATPRODUCT_AtB,C);
9868:   return(0);
9869: }

9871: /*@
9872:    MatMatMatMult - Performs Matrix-Matrix-Matrix Multiplication D=A*B*C.

9874:    Neighbor-wise Collective on Mat

9876:    Input Parameters:
9877: +  A - the left matrix
9878: .  B - the middle matrix
9879: .  C - the right matrix
9880: .  scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
9881: -  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
9882:           if the result is a dense matrix this is irrelevent

9884:    Output Parameters:
9885: .  D - the product matrix

9887:    Notes:
9888:    Unless scall is MAT_REUSE_MATRIX D will be created.

9890:    MAT_REUSE_MATRIX can only be used if the matrices A, B and C have the same nonzero pattern as in the previous call

9892:    To determine the correct fill value, run with -info and search for the string "Fill ratio" to see the value
9893:    actually needed.

9895:    If you have many matrices with the same non-zero structure to multiply, you
9896:    should use MAT_REUSE_MATRIX in all calls but the first or

9898:    Level: intermediate

9900: .seealso: MatMatMult, MatPtAP()
9901: @*/
9902: PetscErrorCode MatMatMatMult(Mat A,Mat B,Mat C,MatReuse scall,PetscReal fill,Mat *D)
9903: {

9907:   if (scall == MAT_REUSE_MATRIX) MatCheckProduct(*D,6);
9908:   if (scall == MAT_INPLACE_MATRIX) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Inplace product not supported");

9910:   if (scall == MAT_INITIAL_MATRIX) {
9911:     MatProductCreate(A,B,C,D);
9912:     MatProductSetType(*D,MATPRODUCT_ABC);
9913:     MatProductSetAlgorithm(*D,"default");
9914:     MatProductSetFill(*D,fill);

9916:     (*D)->product->api_user = PETSC_TRUE;
9917:     MatProductSetFromOptions(*D);
9918:     if (!(*D)->ops->productsymbolic) SETERRQ4(PetscObjectComm((PetscObject)(*D)),PETSC_ERR_SUP,"MatProduct %s not supported for A %s, B %s and C %s",MatProductTypes[MATPRODUCT_ABC],((PetscObject)A)->type_name,((PetscObject)B)->type_name,((PetscObject)C)->type_name);
9919:     MatProductSymbolic(*D);
9920:   } else { /* user may change input matrices when REUSE */
9921:     MatProductReplaceMats(A,B,C,*D);
9922:   }
9923:   MatProductNumeric(*D);
9924:   return(0);
9925: }

9927: /*@
9928:    MatCreateRedundantMatrix - Create redundant matrices and put them into processors of subcommunicators.

9930:    Collective on Mat

9932:    Input Parameters:
9933: +  mat - the matrix
9934: .  nsubcomm - the number of subcommunicators (= number of redundant parallel or sequential matrices)
9935: .  subcomm - MPI communicator split from the communicator where mat resides in (or MPI_COMM_NULL if nsubcomm is used)
9936: -  reuse - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX

9938:    Output Parameter:
9939: .  matredundant - redundant matrix

9941:    Notes:
9942:    MAT_REUSE_MATRIX can only be used when the nonzero structure of the
9943:    original matrix has not changed from that last call to MatCreateRedundantMatrix().

9945:    This routine creates the duplicated matrices in subcommunicators; you should NOT create them before
9946:    calling it.

9948:    Level: advanced


9951: .seealso: MatDestroy()
9952: @*/
9953: PetscErrorCode MatCreateRedundantMatrix(Mat mat,PetscInt nsubcomm,MPI_Comm subcomm,MatReuse reuse,Mat *matredundant)
9954: {
9956:   MPI_Comm       comm;
9957:   PetscMPIInt    size;
9958:   PetscInt       mloc_sub,nloc_sub,rstart,rend,M=mat->rmap->N,N=mat->cmap->N,bs=mat->rmap->bs;
9959:   Mat_Redundant  *redund=NULL;
9960:   PetscSubcomm   psubcomm=NULL;
9961:   MPI_Comm       subcomm_in=subcomm;
9962:   Mat            *matseq;
9963:   IS             isrow,iscol;
9964:   PetscBool      newsubcomm=PETSC_FALSE;

9968:   if (nsubcomm && reuse == MAT_REUSE_MATRIX) {
9971:   }

9973:   MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);
9974:   if (size == 1 || nsubcomm == 1) {
9975:     if (reuse == MAT_INITIAL_MATRIX) {
9976:       MatDuplicate(mat,MAT_COPY_VALUES,matredundant);
9977:     } else {
9978:       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");
9979:       MatCopy(mat,*matredundant,SAME_NONZERO_PATTERN);
9980:     }
9981:     return(0);
9982:   }

9984:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
9985:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
9986:   MatCheckPreallocated(mat,1);

9988:   PetscLogEventBegin(MAT_RedundantMat,mat,0,0,0);
9989:   if (subcomm_in == MPI_COMM_NULL && reuse == MAT_INITIAL_MATRIX) { /* get subcomm if user does not provide subcomm */
9990:     /* create psubcomm, then get subcomm */
9991:     PetscObjectGetComm((PetscObject)mat,&comm);
9992:     MPI_Comm_size(comm,&size);
9993:     if (nsubcomm < 1 || nsubcomm > size) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"nsubcomm must between 1 and %D",size);

9995:     PetscSubcommCreate(comm,&psubcomm);
9996:     PetscSubcommSetNumber(psubcomm,nsubcomm);
9997:     PetscSubcommSetType(psubcomm,PETSC_SUBCOMM_CONTIGUOUS);
9998:     PetscSubcommSetFromOptions(psubcomm);
9999:     PetscCommDuplicate(PetscSubcommChild(psubcomm),&subcomm,NULL);
10000:     newsubcomm = PETSC_TRUE;
10001:     PetscSubcommDestroy(&psubcomm);
10002:   }

10004:   /* get isrow, iscol and a local sequential matrix matseq[0] */
10005:   if (reuse == MAT_INITIAL_MATRIX) {
10006:     mloc_sub = PETSC_DECIDE;
10007:     nloc_sub = PETSC_DECIDE;
10008:     if (bs < 1) {
10009:       PetscSplitOwnership(subcomm,&mloc_sub,&M);
10010:       PetscSplitOwnership(subcomm,&nloc_sub,&N);
10011:     } else {
10012:       PetscSplitOwnershipBlock(subcomm,bs,&mloc_sub,&M);
10013:       PetscSplitOwnershipBlock(subcomm,bs,&nloc_sub,&N);
10014:     }
10015:     MPI_Scan(&mloc_sub,&rend,1,MPIU_INT,MPI_SUM,subcomm);
10016:     rstart = rend - mloc_sub;
10017:     ISCreateStride(PETSC_COMM_SELF,mloc_sub,rstart,1,&isrow);
10018:     ISCreateStride(PETSC_COMM_SELF,N,0,1,&iscol);
10019:   } else { /* reuse == MAT_REUSE_MATRIX */
10020:     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");
10021:     /* retrieve subcomm */
10022:     PetscObjectGetComm((PetscObject)(*matredundant),&subcomm);
10023:     redund = (*matredundant)->redundant;
10024:     isrow  = redund->isrow;
10025:     iscol  = redund->iscol;
10026:     matseq = redund->matseq;
10027:   }
10028:   MatCreateSubMatrices(mat,1,&isrow,&iscol,reuse,&matseq);

10030:   /* get matredundant over subcomm */
10031:   if (reuse == MAT_INITIAL_MATRIX) {
10032:     MatCreateMPIMatConcatenateSeqMat(subcomm,matseq[0],nloc_sub,reuse,matredundant);

10034:     /* create a supporting struct and attach it to C for reuse */
10035:     PetscNewLog(*matredundant,&redund);
10036:     (*matredundant)->redundant = redund;
10037:     redund->isrow              = isrow;
10038:     redund->iscol              = iscol;
10039:     redund->matseq             = matseq;
10040:     if (newsubcomm) {
10041:       redund->subcomm          = subcomm;
10042:     } else {
10043:       redund->subcomm          = MPI_COMM_NULL;
10044:     }
10045:   } else {
10046:     MatCreateMPIMatConcatenateSeqMat(subcomm,matseq[0],PETSC_DECIDE,reuse,matredundant);
10047:   }
10048:   PetscLogEventEnd(MAT_RedundantMat,mat,0,0,0);
10049:   return(0);
10050: }

10052: /*@C
10053:    MatGetMultiProcBlock - Create multiple [bjacobi] 'parallel submatrices' from
10054:    a given 'mat' object. Each submatrix can span multiple procs.

10056:    Collective on Mat

10058:    Input Parameters:
10059: +  mat - the matrix
10060: .  subcomm - the subcommunicator obtained by com_split(comm)
10061: -  scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX

10063:    Output Parameter:
10064: .  subMat - 'parallel submatrices each spans a given subcomm

10066:   Notes:
10067:   The submatrix partition across processors is dictated by 'subComm' a
10068:   communicator obtained by com_split(comm). The comm_split
10069:   is not restriced to be grouped with consecutive original ranks.

10071:   Due the comm_split() usage, the parallel layout of the submatrices
10072:   map directly to the layout of the original matrix [wrt the local
10073:   row,col partitioning]. So the original 'DiagonalMat' naturally maps
10074:   into the 'DiagonalMat' of the subMat, hence it is used directly from
10075:   the subMat. However the offDiagMat looses some columns - and this is
10076:   reconstructed with MatSetValues()

10078:   Level: advanced


10081: .seealso: MatCreateSubMatrices()
10082: @*/
10083: PetscErrorCode   MatGetMultiProcBlock(Mat mat, MPI_Comm subComm, MatReuse scall,Mat *subMat)
10084: {
10086:   PetscMPIInt    commsize,subCommSize;

10089:   MPI_Comm_size(PetscObjectComm((PetscObject)mat),&commsize);
10090:   MPI_Comm_size(subComm,&subCommSize);
10091:   if (subCommSize > commsize) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_OUTOFRANGE,"CommSize %D < SubCommZize %D",commsize,subCommSize);

10093:   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");
10094:   PetscLogEventBegin(MAT_GetMultiProcBlock,mat,0,0,0);
10095:   (*mat->ops->getmultiprocblock)(mat,subComm,scall,subMat);
10096:   PetscLogEventEnd(MAT_GetMultiProcBlock,mat,0,0,0);
10097:   return(0);
10098: }

10100: /*@
10101:    MatGetLocalSubMatrix - Gets a reference to a submatrix specified in local numbering

10103:    Not Collective

10105:    Input Arguments:
10106: +  mat - matrix to extract local submatrix from
10107: .  isrow - local row indices for submatrix
10108: -  iscol - local column indices for submatrix

10110:    Output Arguments:
10111: .  submat - the submatrix

10113:    Level: intermediate

10115:    Notes:
10116:    The submat should be returned with MatRestoreLocalSubMatrix().

10118:    Depending on the format of mat, the returned submat may not implement MatMult().  Its communicator may be
10119:    the same as mat, it may be PETSC_COMM_SELF, or some other subcomm of mat's.

10121:    The submat always implements MatSetValuesLocal().  If isrow and iscol have the same block size, then
10122:    MatSetValuesBlockedLocal() will also be implemented.

10124:    The mat must have had a ISLocalToGlobalMapping provided to it with MatSetLocalToGlobalMapping(). Note that
10125:    matrices obtained with DMCreateMatrix() generally already have the local to global mapping provided.

10127: .seealso: MatRestoreLocalSubMatrix(), MatCreateLocalRef(), MatSetLocalToGlobalMapping()
10128: @*/
10129: PetscErrorCode MatGetLocalSubMatrix(Mat mat,IS isrow,IS iscol,Mat *submat)
10130: {

10139:   if (!mat->rmap->mapping) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Matrix must have local to global mapping provided before this call");

10141:   if (mat->ops->getlocalsubmatrix) {
10142:     (*mat->ops->getlocalsubmatrix)(mat,isrow,iscol,submat);
10143:   } else {
10144:     MatCreateLocalRef(mat,isrow,iscol,submat);
10145:   }
10146:   return(0);
10147: }

10149: /*@
10150:    MatRestoreLocalSubMatrix - Restores a reference to a submatrix specified in local numbering

10152:    Not Collective

10154:    Input Arguments:
10155:    mat - matrix to extract local submatrix from
10156:    isrow - local row indices for submatrix
10157:    iscol - local column indices for submatrix
10158:    submat - the submatrix

10160:    Level: intermediate

10162: .seealso: MatGetLocalSubMatrix()
10163: @*/
10164: PetscErrorCode MatRestoreLocalSubMatrix(Mat mat,IS isrow,IS iscol,Mat *submat)
10165: {

10174:   if (*submat) {
10176:   }

10178:   if (mat->ops->restorelocalsubmatrix) {
10179:     (*mat->ops->restorelocalsubmatrix)(mat,isrow,iscol,submat);
10180:   } else {
10181:     MatDestroy(submat);
10182:   }
10183:   *submat = NULL;
10184:   return(0);
10185: }

10187: /* --------------------------------------------------------*/
10188: /*@
10189:    MatFindZeroDiagonals - Finds all the rows of a matrix that have zero or no diagonal entry in the matrix

10191:    Collective on Mat

10193:    Input Parameter:
10194: .  mat - the matrix

10196:    Output Parameter:
10197: .  is - if any rows have zero diagonals this contains the list of them

10199:    Level: developer

10201: .seealso: MatMultTranspose(), MatMultAdd(), MatMultTransposeAdd()
10202: @*/
10203: PetscErrorCode MatFindZeroDiagonals(Mat mat,IS *is)
10204: {

10210:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
10211:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");

10213:   if (!mat->ops->findzerodiagonals) {
10214:     Vec                diag;
10215:     const PetscScalar *a;
10216:     PetscInt          *rows;
10217:     PetscInt           rStart, rEnd, r, nrow = 0;

10219:     MatCreateVecs(mat, &diag, NULL);
10220:     MatGetDiagonal(mat, diag);
10221:     MatGetOwnershipRange(mat, &rStart, &rEnd);
10222:     VecGetArrayRead(diag, &a);
10223:     for (r = 0; r < rEnd-rStart; ++r) if (a[r] == 0.0) ++nrow;
10224:     PetscMalloc1(nrow, &rows);
10225:     nrow = 0;
10226:     for (r = 0; r < rEnd-rStart; ++r) if (a[r] == 0.0) rows[nrow++] = r+rStart;
10227:     VecRestoreArrayRead(diag, &a);
10228:     VecDestroy(&diag);
10229:     ISCreateGeneral(PetscObjectComm((PetscObject) mat), nrow, rows, PETSC_OWN_POINTER, is);
10230:   } else {
10231:     (*mat->ops->findzerodiagonals)(mat, is);
10232:   }
10233:   return(0);
10234: }

10236: /*@
10237:    MatFindOffBlockDiagonalEntries - Finds all the rows of a matrix that have entries outside of the main diagonal block (defined by the matrix block size)

10239:    Collective on Mat

10241:    Input Parameter:
10242: .  mat - the matrix

10244:    Output Parameter:
10245: .  is - contains the list of rows with off block diagonal entries

10247:    Level: developer

10249: .seealso: MatMultTranspose(), MatMultAdd(), MatMultTransposeAdd()
10250: @*/
10251: PetscErrorCode MatFindOffBlockDiagonalEntries(Mat mat,IS *is)
10252: {

10258:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
10259:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");

10261:   if (!mat->ops->findoffblockdiagonalentries) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Matrix type %s does not have a find off block diagonal entries defined",((PetscObject)mat)->type_name);
10262:   (*mat->ops->findoffblockdiagonalentries)(mat,is);
10263:   return(0);
10264: }

10266: /*@C
10267:   MatInvertBlockDiagonal - Inverts the block diagonal entries.

10269:   Collective on Mat

10271:   Input Parameters:
10272: . mat - the matrix

10274:   Output Parameters:
10275: . values - the block inverses in column major order (FORTRAN-like)

10277:    Note:
10278:    This routine is not available from Fortran.

10280:   Level: advanced

10282: .seealso: MatInvertBockDiagonalMat
10283: @*/
10284: PetscErrorCode MatInvertBlockDiagonal(Mat mat,const PetscScalar **values)
10285: {

10290:   if (!mat->assembled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
10291:   if (mat->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
10292:   if (!mat->ops->invertblockdiagonal) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not supported for type %s",((PetscObject)mat)->type_name);
10293:   (*mat->ops->invertblockdiagonal)(mat,values);
10294:   return(0);
10295: }

10297: /*@C
10298:   MatInvertVariableBlockDiagonal - Inverts the block diagonal entries.

10300:   Collective on Mat

10302:   Input Parameters:
10303: + mat - the matrix
10304: . nblocks - the number of blocks
10305: - bsizes - the size of each block

10307:   Output Parameters:
10308: . values - the block inverses in column major order (FORTRAN-like)

10310:    Note:
10311:    This routine is not available from Fortran.

10313:   Level: advanced

10315: .seealso: MatInvertBockDiagonal()
10316: @*/
10317: PetscErrorCode MatInvertVariableBlockDiagonal(Mat mat,PetscInt nblocks,const PetscInt *bsizes,PetscScalar *values)
10318: {

10323:   if (!mat->assembled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
10324:   if (mat->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
10325:   if (!mat->ops->invertvariableblockdiagonal) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not supported for type",((PetscObject)mat)->type_name);
10326:   (*mat->ops->invertvariableblockdiagonal)(mat,nblocks,bsizes,values);
10327:   return(0);
10328: }

10330: /*@
10331:   MatInvertBlockDiagonalMat - set matrix C to be the inverted block diagonal of matrix A

10333:   Collective on Mat

10335:   Input Parameters:
10336: . A - the matrix

10338:   Output Parameters:
10339: . C - matrix with inverted block diagonal of A.  This matrix should be created and may have its type set.

10341:   Notes: the blocksize of the matrix is used to determine the blocks on the diagonal of C

10343:   Level: advanced

10345: .seealso: MatInvertBockDiagonal()
10346: @*/
10347: PetscErrorCode MatInvertBlockDiagonalMat(Mat A,Mat C)
10348: {
10349:   PetscErrorCode     ierr;
10350:   const PetscScalar *vals;
10351:   PetscInt          *dnnz;
10352:   PetscInt           M,N,m,n,rstart,rend,bs,i,j;

10355:   MatInvertBlockDiagonal(A,&vals);
10356:   MatGetBlockSize(A,&bs);
10357:   MatGetSize(A,&M,&N);
10358:   MatGetLocalSize(A,&m,&n);
10359:   MatSetSizes(C,m,n,M,N);
10360:   MatSetBlockSize(C,bs);
10361:   PetscMalloc1(m/bs,&dnnz);
10362:   for (j = 0; j < m/bs; j++) dnnz[j] = 1;
10363:   MatXAIJSetPreallocation(C,bs,dnnz,NULL,NULL,NULL);
10364:   PetscFree(dnnz);
10365:   MatGetOwnershipRange(C,&rstart,&rend);
10366:   MatSetOption(C,MAT_ROW_ORIENTED,PETSC_FALSE);
10367:   for (i = rstart/bs; i < rend/bs; i++) {
10368:     MatSetValuesBlocked(C,1,&i,1,&i,&vals[(i-rstart/bs)*bs*bs],INSERT_VALUES);
10369:   }
10370:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
10371:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
10372:   MatSetOption(C,MAT_ROW_ORIENTED,PETSC_TRUE);
10373:   return(0);
10374: }

10376: /*@C
10377:     MatTransposeColoringDestroy - Destroys a coloring context for matrix product C=A*B^T that was created
10378:     via MatTransposeColoringCreate().

10380:     Collective on MatTransposeColoring

10382:     Input Parameter:
10383: .   c - coloring context

10385:     Level: intermediate

10387: .seealso: MatTransposeColoringCreate()
10388: @*/
10389: PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring *c)
10390: {
10391:   PetscErrorCode       ierr;
10392:   MatTransposeColoring matcolor=*c;

10395:   if (!matcolor) return(0);
10396:   if (--((PetscObject)matcolor)->refct > 0) {matcolor = NULL; return(0);}

10398:   PetscFree3(matcolor->ncolumns,matcolor->nrows,matcolor->colorforrow);
10399:   PetscFree(matcolor->rows);
10400:   PetscFree(matcolor->den2sp);
10401:   PetscFree(matcolor->colorforcol);
10402:   PetscFree(matcolor->columns);
10403:   if (matcolor->brows>0) {
10404:     PetscFree(matcolor->lstart);
10405:   }
10406:   PetscHeaderDestroy(c);
10407:   return(0);
10408: }

10410: /*@C
10411:     MatTransColoringApplySpToDen - Given a symbolic matrix product C=A*B^T for which
10412:     a MatTransposeColoring context has been created, computes a dense B^T by Apply
10413:     MatTransposeColoring to sparse B.

10415:     Collective on MatTransposeColoring

10417:     Input Parameters:
10418: +   B - sparse matrix B
10419: .   Btdense - symbolic dense matrix B^T
10420: -   coloring - coloring context created with MatTransposeColoringCreate()

10422:     Output Parameter:
10423: .   Btdense - dense matrix B^T

10425:     Level: advanced

10427:      Notes:
10428:     These are used internally for some implementations of MatRARt()

10430: .seealso: MatTransposeColoringCreate(), MatTransposeColoringDestroy(), MatTransColoringApplyDenToSp()

10432: @*/
10433: PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring coloring,Mat B,Mat Btdense)
10434: {


10442:   if (!B->ops->transcoloringapplysptoden) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not supported for this matrix type %s",((PetscObject)B)->type_name);
10443:   (B->ops->transcoloringapplysptoden)(coloring,B,Btdense);
10444:   return(0);
10445: }

10447: /*@C
10448:     MatTransColoringApplyDenToSp - Given a symbolic matrix product Csp=A*B^T for which
10449:     a MatTransposeColoring context has been created and a dense matrix Cden=A*Btdense
10450:     in which Btdens is obtained from MatTransColoringApplySpToDen(), recover sparse matrix
10451:     Csp from Cden.

10453:     Collective on MatTransposeColoring

10455:     Input Parameters:
10456: +   coloring - coloring context created with MatTransposeColoringCreate()
10457: -   Cden - matrix product of a sparse matrix and a dense matrix Btdense

10459:     Output Parameter:
10460: .   Csp - sparse matrix

10462:     Level: advanced

10464:      Notes:
10465:     These are used internally for some implementations of MatRARt()

10467: .seealso: MatTransposeColoringCreate(), MatTransposeColoringDestroy(), MatTransColoringApplySpToDen()

10469: @*/
10470: PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring matcoloring,Mat Cden,Mat Csp)
10471: {


10479:   if (!Csp->ops->transcoloringapplydentosp) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not supported for this matrix type %s",((PetscObject)Csp)->type_name);
10480:   (Csp->ops->transcoloringapplydentosp)(matcoloring,Cden,Csp);
10481:   MatAssemblyBegin(Csp,MAT_FINAL_ASSEMBLY);
10482:   MatAssemblyEnd(Csp,MAT_FINAL_ASSEMBLY);
10483:   return(0);
10484: }

10486: /*@C
10487:    MatTransposeColoringCreate - Creates a matrix coloring context for matrix product C=A*B^T.

10489:    Collective on Mat

10491:    Input Parameters:
10492: +  mat - the matrix product C
10493: -  iscoloring - the coloring of the matrix; usually obtained with MatColoringCreate() or DMCreateColoring()

10495:     Output Parameter:
10496: .   color - the new coloring context

10498:     Level: intermediate

10500: .seealso: MatTransposeColoringDestroy(),  MatTransColoringApplySpToDen(),
10501:            MatTransColoringApplyDenToSp()
10502: @*/
10503: PetscErrorCode MatTransposeColoringCreate(Mat mat,ISColoring iscoloring,MatTransposeColoring *color)
10504: {
10505:   MatTransposeColoring c;
10506:   MPI_Comm             comm;
10507:   PetscErrorCode       ierr;

10510:   PetscLogEventBegin(MAT_TransposeColoringCreate,mat,0,0,0);
10511:   PetscObjectGetComm((PetscObject)mat,&comm);
10512:   PetscHeaderCreate(c,MAT_TRANSPOSECOLORING_CLASSID,"MatTransposeColoring","Matrix product C=A*B^T via coloring","Mat",comm,MatTransposeColoringDestroy,NULL);

10514:   c->ctype = iscoloring->ctype;
10515:   if (mat->ops->transposecoloringcreate) {
10516:     (*mat->ops->transposecoloringcreate)(mat,iscoloring,c);
10517:   } else SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Code not yet written for matrix type %s",((PetscObject)mat)->type_name);

10519:   *color = c;
10520:   PetscLogEventEnd(MAT_TransposeColoringCreate,mat,0,0,0);
10521:   return(0);
10522: }

10524: /*@
10525:       MatGetNonzeroState - Returns a 64 bit integer representing the current state of nonzeros in the matrix. If the
10526:         matrix has had no new nonzero locations added to the matrix since the previous call then the value will be the
10527:         same, otherwise it will be larger

10529:      Not Collective

10531:   Input Parameter:
10532: .    A  - the matrix

10534:   Output Parameter:
10535: .    state - the current state

10537:   Notes:
10538:     You can only compare states from two different calls to the SAME matrix, you cannot compare calls between
10539:          different matrices

10541:   Level: intermediate

10543: @*/
10544: PetscErrorCode MatGetNonzeroState(Mat mat,PetscObjectState *state)
10545: {
10548:   *state = mat->nonzerostate;
10549:   return(0);
10550: }

10552: /*@
10553:       MatCreateMPIMatConcatenateSeqMat - Creates a single large PETSc matrix by concatenating sequential
10554:                  matrices from each processor

10556:     Collective

10558:    Input Parameters:
10559: +    comm - the communicators the parallel matrix will live on
10560: .    seqmat - the input sequential matrices
10561: .    n - number of local columns (or PETSC_DECIDE)
10562: -    reuse - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX

10564:    Output Parameter:
10565: .    mpimat - the parallel matrix generated

10567:     Level: advanced

10569:    Notes:
10570:     The number of columns of the matrix in EACH processor MUST be the same.

10572: @*/
10573: PetscErrorCode MatCreateMPIMatConcatenateSeqMat(MPI_Comm comm,Mat seqmat,PetscInt n,MatReuse reuse,Mat *mpimat)
10574: {

10578:   if (!seqmat->ops->creatempimatconcatenateseqmat) SETERRQ1(PetscObjectComm((PetscObject)seqmat),PETSC_ERR_SUP,"Mat type %s",((PetscObject)seqmat)->type_name);
10579:   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");

10581:   PetscLogEventBegin(MAT_Merge,seqmat,0,0,0);
10582:   (*seqmat->ops->creatempimatconcatenateseqmat)(comm,seqmat,n,reuse,mpimat);
10583:   PetscLogEventEnd(MAT_Merge,seqmat,0,0,0);
10584:   return(0);
10585: }

10587: /*@
10588:      MatSubdomainsCreateCoalesce - Creates index subdomains by coalescing adjacent
10589:                  ranks' ownership ranges.

10591:     Collective on A

10593:    Input Parameters:
10594: +    A   - the matrix to create subdomains from
10595: -    N   - requested number of subdomains


10598:    Output Parameters:
10599: +    n   - number of subdomains resulting on this rank
10600: -    iss - IS list with indices of subdomains on this rank

10602:     Level: advanced

10604:     Notes:
10605:     number of subdomains must be smaller than the communicator size
10606: @*/
10607: PetscErrorCode MatSubdomainsCreateCoalesce(Mat A,PetscInt N,PetscInt *n,IS *iss[])
10608: {
10609:   MPI_Comm        comm,subcomm;
10610:   PetscMPIInt     size,rank,color;
10611:   PetscInt        rstart,rend,k;
10612:   PetscErrorCode  ierr;

10615:   PetscObjectGetComm((PetscObject)A,&comm);
10616:   MPI_Comm_size(comm,&size);
10617:   MPI_Comm_rank(comm,&rank);
10618:   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);
10619:   *n = 1;
10620:   k = ((PetscInt)size)/N + ((PetscInt)size%N>0); /* There are up to k ranks to a color */
10621:   color = rank/k;
10622:   MPI_Comm_split(comm,color,rank,&subcomm);
10623:   PetscMalloc1(1,iss);
10624:   MatGetOwnershipRange(A,&rstart,&rend);
10625:   ISCreateStride(subcomm,rend-rstart,rstart,1,iss[0]);
10626:   MPI_Comm_free(&subcomm);
10627:   return(0);
10628: }

10630: /*@
10631:    MatGalerkin - Constructs the coarse grid problem via Galerkin projection.

10633:    If the interpolation and restriction operators are the same, uses MatPtAP.
10634:    If they are not the same, use MatMatMatMult.

10636:    Once the coarse grid problem is constructed, correct for interpolation operators
10637:    that are not of full rank, which can legitimately happen in the case of non-nested
10638:    geometric multigrid.

10640:    Input Parameters:
10641: +  restrct - restriction operator
10642: .  dA - fine grid matrix
10643: .  interpolate - interpolation operator
10644: .  reuse - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
10645: -  fill - expected fill, use PETSC_DEFAULT if you do not have a good estimate

10647:    Output Parameters:
10648: .  A - the Galerkin coarse matrix

10650:    Options Database Key:
10651: .  -pc_mg_galerkin <both,pmat,mat,none>

10653:    Level: developer

10655: .seealso: MatPtAP(), MatMatMatMult()
10656: @*/
10657: PetscErrorCode  MatGalerkin(Mat restrct, Mat dA, Mat interpolate, MatReuse reuse, PetscReal fill, Mat *A)
10658: {
10660:   IS             zerorows;
10661:   Vec            diag;

10664:   if (reuse == MAT_INPLACE_MATRIX) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Inplace product not supported");
10665:   /* Construct the coarse grid matrix */
10666:   if (interpolate == restrct) {
10667:     MatPtAP(dA,interpolate,reuse,fill,A);
10668:   } else {
10669:     MatMatMatMult(restrct,dA,interpolate,reuse,fill,A);
10670:   }

10672:   /* If the interpolation matrix is not of full rank, A will have zero rows.
10673:      This can legitimately happen in the case of non-nested geometric multigrid.
10674:      In that event, we set the rows of the matrix to the rows of the identity,
10675:      ignoring the equations (as the RHS will also be zero). */

10677:   MatFindZeroRows(*A, &zerorows);

10679:   if (zerorows != NULL) { /* if there are any zero rows */
10680:     MatCreateVecs(*A, &diag, NULL);
10681:     MatGetDiagonal(*A, diag);
10682:     VecISSet(diag, zerorows, 1.0);
10683:     MatDiagonalSet(*A, diag, INSERT_VALUES);
10684:     VecDestroy(&diag);
10685:     ISDestroy(&zerorows);
10686:   }
10687:   return(0);
10688: }

10690: /*@C
10691:     MatSetOperation - Allows user to set a matrix operation for any matrix type

10693:    Logically Collective on Mat

10695:     Input Parameters:
10696: +   mat - the matrix
10697: .   op - the name of the operation
10698: -   f - the function that provides the operation

10700:    Level: developer

10702:     Usage:
10703: $      extern PetscErrorCode usermult(Mat,Vec,Vec);
10704: $      MatCreateXXX(comm,...&A);
10705: $      MatSetOperation(A,MATOP_MULT,(void(*)(void))usermult);

10707:     Notes:
10708:     See the file include/petscmat.h for a complete list of matrix
10709:     operations, which all have the form MATOP_<OPERATION>, where
10710:     <OPERATION> is the name (in all capital letters) of the
10711:     user interface routine (e.g., MatMult() -> MATOP_MULT).

10713:     All user-provided functions (except for MATOP_DESTROY) should have the same calling
10714:     sequence as the usual matrix interface routines, since they
10715:     are intended to be accessed via the usual matrix interface
10716:     routines, e.g.,
10717: $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)

10719:     In particular each function MUST return an error code of 0 on success and
10720:     nonzero on failure.

10722:     This routine is distinct from MatShellSetOperation() in that it can be called on any matrix type.

10724: .seealso: MatGetOperation(), MatCreateShell(), MatShellSetContext(), MatShellSetOperation()
10725: @*/
10726: PetscErrorCode MatSetOperation(Mat mat,MatOperation op,void (*f)(void))
10727: {
10730:   if (op == MATOP_VIEW && !mat->ops->viewnative && f != (void (*)(void))(mat->ops->view)) {
10731:     mat->ops->viewnative = mat->ops->view;
10732:   }
10733:   (((void(**)(void))mat->ops)[op]) = f;
10734:   return(0);
10735: }

10737: /*@C
10738:     MatGetOperation - Gets a matrix operation for any matrix type.

10740:     Not Collective

10742:     Input Parameters:
10743: +   mat - the matrix
10744: -   op - the name of the operation

10746:     Output Parameter:
10747: .   f - the function that provides the operation

10749:     Level: developer

10751:     Usage:
10752: $      PetscErrorCode (*usermult)(Mat,Vec,Vec);
10753: $      MatGetOperation(A,MATOP_MULT,(void(**)(void))&usermult);

10755:     Notes:
10756:     See the file include/petscmat.h for a complete list of matrix
10757:     operations, which all have the form MATOP_<OPERATION>, where
10758:     <OPERATION> is the name (in all capital letters) of the
10759:     user interface routine (e.g., MatMult() -> MATOP_MULT).

10761:     This routine is distinct from MatShellGetOperation() in that it can be called on any matrix type.

10763: .seealso: MatSetOperation(), MatCreateShell(), MatShellGetContext(), MatShellGetOperation()
10764: @*/
10765: PetscErrorCode MatGetOperation(Mat mat,MatOperation op,void(**f)(void))
10766: {
10769:   *f = (((void (**)(void))mat->ops)[op]);
10770:   return(0);
10771: }

10773: /*@
10774:     MatHasOperation - Determines whether the given matrix supports the particular
10775:     operation.

10777:    Not Collective

10779:    Input Parameters:
10780: +  mat - the matrix
10781: -  op - the operation, for example, MATOP_GET_DIAGONAL

10783:    Output Parameter:
10784: .  has - either PETSC_TRUE or PETSC_FALSE

10786:    Level: advanced

10788:    Notes:
10789:    See the file include/petscmat.h for a complete list of matrix
10790:    operations, which all have the form MATOP_<OPERATION>, where
10791:    <OPERATION> is the name (in all capital letters) of the
10792:    user-level routine.  E.g., MatNorm() -> MATOP_NORM.

10794: .seealso: MatCreateShell()
10795: @*/
10796: PetscErrorCode MatHasOperation(Mat mat,MatOperation op,PetscBool *has)
10797: {

10802:   /* symbolic product can be set before matrix type */
10805:   if (mat->ops->hasoperation) {
10806:     (*mat->ops->hasoperation)(mat,op,has);
10807:   } else {
10808:     if (((void**)mat->ops)[op]) *has =  PETSC_TRUE;
10809:     else {
10810:       *has = PETSC_FALSE;
10811:       if (op == MATOP_CREATE_SUBMATRIX) {
10812:         PetscMPIInt size;

10814:         MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);
10815:         if (size == 1) {
10816:           MatHasOperation(mat,MATOP_CREATE_SUBMATRICES,has);
10817:         }
10818:       }
10819:     }
10820:   }
10821:   return(0);
10822: }

10824: /*@
10825:     MatHasCongruentLayouts - Determines whether the rows and columns layouts
10826:     of the matrix are congruent

10828:    Collective on mat

10830:    Input Parameters:
10831: .  mat - the matrix

10833:    Output Parameter:
10834: .  cong - either PETSC_TRUE or PETSC_FALSE

10836:    Level: beginner

10838:    Notes:

10840: .seealso: MatCreate(), MatSetSizes()
10841: @*/
10842: PetscErrorCode MatHasCongruentLayouts(Mat mat,PetscBool *cong)
10843: {

10850:   if (!mat->rmap || !mat->cmap) {
10851:     *cong = mat->rmap == mat->cmap ? PETSC_TRUE : PETSC_FALSE;
10852:     return(0);
10853:   }
10854:   if (mat->congruentlayouts == PETSC_DECIDE) { /* first time we compare rows and cols layouts */
10855:     PetscLayoutCompare(mat->rmap,mat->cmap,cong);
10856:     if (*cong) mat->congruentlayouts = 1;
10857:     else       mat->congruentlayouts = 0;
10858:   } else *cong = mat->congruentlayouts ? PETSC_TRUE : PETSC_FALSE;
10859:   return(0);
10860: }

10862: PetscErrorCode MatSetInf(Mat A)
10863: {

10867:   if (!A->ops->setinf) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for this operation for this matrix type");
10868:   (*A->ops->setinf)(A);
10869:   return(0);
10870: }