xref: /petsc/src/mat/interface/matrix.c (revision 7a2efc394589b136789e5cec3f47c1a804ea8ac0)
1 
2 #ifndef lint
3 static char vcid[] = "$Id: matrix.c,v 1.90 1995/10/01 21:52:25 bsmith Exp curfman $";
4 #endif
5 
6 /*
7    This is where the abstract matrix operations are defined
8 */
9 
10 #include "petsc.h"
11 #include "matimpl.h"        /*I "mat.h" I*/
12 #include "vec/vecimpl.h"
13 #include "pinclude/pviewer.h"
14 
15 extern int MatGetReorderingMethodFromOptions_Private(MatOrdering *);
16 
17 /*@C
18    MatGetReordering - Gets a reordering for a matrix to reduce fill or to
19    improve numerical stability of LU factorization.
20 
21    Input Parameters:
22 .  mat - the matrix
23 .  type - type of reordering, one of the following:
24 $      ORDER_NATURAL - Natural
25 $      ORDER_ND - Nested Dissection
26 $      ORDER_1WD - One-way Dissection
27 $      ORDER_RCM - Reverse Cuthill-McGee
28 $      ORDER_QMD - Quotient Minimum Degree
29 
30    Output Parameters:
31 .  rperm - row permutation indices
32 .  cperm - column permutation indices
33 
34    Options Database Keys:
35    To specify the ordering through the options database, use one of
36    the following
37 $    -mat_order natural, -mat_order nd, -mat_order 1wd,
38 $    -mat_order rcm, -mat_order qmd
39 
40    Notes:
41    If the column permutations and row permutations are the same,
42    then MatGetReordering() returns 0 in cperm.
43 
44 .keywords: matrix, set, ordering, factorization, direct, ILU, LU,
45            fill, reordering, natural, Nested Dissection,
46            One-way Dissection, Cholesky, Reverse Cuthill-McGee,
47            Quotient Minimum Degree
48 @*/
49 int MatGetReordering(Mat mat,MatOrdering type,IS *rperm,IS *cperm)
50 {
51   int         ierr;
52   MatOrdering newtype;
53   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
54   if (!mat->ops.getreordering) {*rperm = 0; *cperm = 0; return 0;}
55   PLogEventBegin(MAT_GetReordering,mat,0,0,0);
56   if (MatGetReorderingMethodFromOptions_Private(&newtype)) type = newtype;
57   ierr = (*mat->ops.getreordering)(mat,type,rperm,cperm); CHKERRQ(ierr);
58   PLogEventEnd(MAT_GetReordering,mat,0,0,0);
59   return 0;
60 }
61 
62 /*@C
63    MatGetRow - Gets a row of a matrix.  You MUST call MatRestoreRow()
64    for each row that you get to ensure that your application does
65    not bleed memory.
66 
67    Input Parameters:
68 .  mat - the matrix
69 .  row - the row to get
70 
71    Output Parameters:
72 .  ncols -  the number of nonzeros in the row
73 .  cols - if nonzero, the column numbers
74 .  vals - if nonzero, the values
75 
76    Notes:
77    This routine is provided for people who need to have direct access
78    to the structure of a matrix.  We hope that we provide enough
79    high-level matrix routines that few users will need it.
80 
81    For better efficiency, set cols and/or vals to zero if you do not
82    wish to extract these quantities.
83 
84 .keywords: matrix, row, get, extract
85 
86 .seealso: MatRestoreRow()
87 @*/
88 int MatGetRow(Mat mat,int row,int *ncols,int **cols,Scalar **vals)
89 {
90   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
91   return (*mat->ops.getrow)(mat,row,ncols,cols,vals);
92 }
93 
94 /*@C
95    MatRestoreRow - Frees any temporary space allocated by MatGetRow().
96 
97    Input Parameters:
98 .  mat - the matrix
99 .  row - the row to get
100 .  ncols, cols - the number of nonzeros and their columns
101 .  vals - if nonzero the column values
102 
103 .keywords: matrix, row, restore
104 
105 .seealso:  MatGetRow()
106 @*/
107 int MatRestoreRow(Mat mat,int row,int *ncols,int **cols,Scalar **vals)
108 {
109   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
110   if (!mat->ops.restorerow) return 0;
111   return (*mat->ops.restorerow)(mat,row,ncols,cols,vals);
112 }
113 /*@
114    MatView - Visualizes a matrix object.
115 
116    Input Parameters:
117 .  mat - the matrix
118 .  ptr - visualization context
119 
120    Notes:
121    The available visualization contexts include
122 $     STDOUT_VIEWER_SELF - standard output (default)
123 $     STDOUT_VIEWER_WORLD - synchronized standard
124 $       output where only the first processor opens
125 $       the file.  All other processors send their
126 $       data to the first processor to print.
127 
128    The user can open alternative vistualization contexts with
129 $    ViewerFileOpenASCII() - output to a specified file
130 $    DrawOpenX() - output nonzero matrix structure to
131 $         an X window display
132 $    ViewerMatlabOpen() - output matrix to Matlab viewer.
133 $         Currently only the sequential dense and AIJ
134 $         matrix types support the Matlab viewer.
135 
136    The user can call ViewerFileSetFormat() to specify the output
137    format.  Available formats include
138 $    FILE_FORMAT_DEFAULT - default, prints matrix contents
139 $    FILE_FORMAT_IMPL - implementation-specific format
140 $      (which is in many cases the same as the default)
141 $    FILE_FORMAT_INFO - basic information about the matrix
142 $      size and structure (not the matrix entries)
143 
144 .keywords: matrix, view, visualize
145 
146 .seealso: ViewerFileSetFormat(), ViewerFileOpenASCII(), DrawOpenX(),
147           ViewerMatlabOpen()
148 @*/
149 int MatView(Mat mat,Viewer ptr)
150 {
151   int format, ierr, rows, cols,nz, nzalloc, mem;
152   FILE *fd;
153   char *cstring;
154   PetscObject  vobj = (PetscObject) ptr;
155 
156   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
157   if (!ptr) { /* so that viewers may be used from debuggers */
158     ptr = STDOUT_VIEWER_SELF; vobj = (PetscObject) ptr;
159   }
160   ierr = ViewerFileGetFormat_Private(ptr,&format); CHKERRQ(ierr);
161   ierr = ViewerFileGetPointer_Private(ptr,&fd); CHKERRQ(ierr);
162   if (vobj->cookie == VIEWER_COOKIE && format == FILE_FORMAT_INFO &&
163     (vobj->type == ASCII_FILE_VIEWER || vobj->type == ASCII_FILES_VIEWER)) {
164     MPIU_fprintf(mat->comm,fd,"Matrix Object:\n");
165     ierr = MatGetName(mat,&cstring); CHKERRQ(ierr);
166     ierr = MatGetSize(mat,&rows,&cols); CHKERRQ(ierr);
167     MPIU_fprintf(mat->comm,fd,
168       "  type=%s, rows=%d, cols=%d\n",cstring,rows,cols);
169     if (mat->ops.getinfo) {
170       ierr = MatGetInfo(mat,MAT_LOCAL,&nz,&nzalloc,&mem); CHKERRQ(ierr);
171       MPIU_fprintf(mat->comm,fd,
172         "  nonzeros=%d, allocated nonzeros=%d\n",nz,nzalloc);
173     }
174   }
175   if (mat->view) {ierr = (*mat->view)((PetscObject)mat,ptr); CHKERRQ(ierr);}
176   return 0;
177 }
178 /*@C
179    MatDestroy - Frees space taken by a matrix.
180 
181    Input Parameter:
182 .  mat - the matrix
183 
184 .keywords: matrix, destroy
185 @*/
186 int MatDestroy(Mat mat)
187 {
188   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
189   return (*mat->destroy)((PetscObject)mat);
190 }
191 /*@
192    MatValidMatrix - Returns 1 if a valid matrix else 0.
193 
194    Input Parameter:
195 .  m - the matrix to check
196 
197 .keywords: matrix, valid
198 @*/
199 int MatValidMatrix(Mat m)
200 {
201   if (!m) return 0;
202   if (m->cookie != MAT_COOKIE) return 0;
203   return 1;
204 }
205 
206 /*@
207    MatSetValues - Inserts or adds a block of values into a matrix.
208    These values may be cached, so MatAssemblyBegin() and MatAssemblyEnd()
209    MUST be called after all calls to MatSetValues() have been completed.
210 
211    Input Parameters:
212 .  mat - the matrix
213 .  v - a logically two-dimensional array of values
214 .  m, indexm - the number of rows and their global indices
215 .  n, indexn - the number of columns and their global indices
216 .  addv - either ADD_VALUES or INSERT_VALUES, where
217 $     ADD_VALUES - adds values to any existing entries
218 $     INSERT_VALUES - replaces existing entries with new values
219 
220    Notes:
221    By default the values, v, are row-oriented and unsorted.
222    See MatSetOptions() for other options.
223 
224    Calls to MatSetValues() with the INSERT_VALUES and ADD_VALUES
225    options cannot be mixed without intervening calls to the assembly
226    routines.
227 
228 .keywords: matrix, insert, add, set, values
229 
230 .seealso: MatSetOptions(), MatAssemblyBegin(), MatAssemblyEnd()
231 @*/
232 int MatSetValues(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v,
233                                                         InsertMode addv)
234 {
235   int ierr;
236   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
237   PLogEventBegin(MAT_SetValues,mat,0,0,0);
238   ierr = (*mat->ops.setvalues)(mat,m,idxm,n,idxn,v,addv);CHKERRQ(ierr);
239   PLogEventEnd(MAT_SetValues,mat,0,0,0);
240   return 0;
241 }
242 
243 /* --------------------------------------------------------*/
244 /*@
245    MatMult - Computes matrix-vector product.
246 
247    Input Parameters:
248 .  mat - the matrix
249 .  x   - the vector to be multilplied
250 
251    Output Parameters:
252 .  y - the result
253 
254 .keywords: matrix, multiply, matrix-vector product
255 
256 .seealso: MatMultTrans(), MatMultAdd(), MatMultTransAdd()
257 @*/
258 int MatMult(Mat mat,Vec x,Vec y)
259 {
260   int ierr;
261   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
262   PETSCVALIDHEADERSPECIFIC(x,VEC_COOKIE);PETSCVALIDHEADERSPECIFIC(y,VEC_COOKIE);
263   PLogEventBegin(MAT_Mult,mat,x,y,0);
264   ierr = (*mat->ops.mult)(mat,x,y); CHKERRQ(ierr);
265   PLogEventEnd(MAT_Mult,mat,x,y,0);
266   return 0;
267 }
268 /*@
269    MatMultTrans - Computes matrix transpose times a vector.
270 
271    Input Parameters:
272 .  mat - the matrix
273 .  x   - the vector to be multilplied
274 
275    Output Parameters:
276 .  y - the result
277 
278 .keywords: matrix, multiply, matrix-vector product, transpose
279 
280 .seealso: MatMult(), MatMultAdd(), MatMultTransAdd()
281 @*/
282 int MatMultTrans(Mat mat,Vec x,Vec y)
283 {
284   int ierr;
285   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
286   PETSCVALIDHEADERSPECIFIC(x,VEC_COOKIE); PETSCVALIDHEADERSPECIFIC(y,VEC_COOKIE);
287   PLogEventBegin(MAT_MultTrans,mat,x,y,0);
288   ierr = (*mat->ops.multtrans)(mat,x,y); CHKERRQ(ierr);
289   PLogEventEnd(MAT_MultTrans,mat,x,y,0);
290   return 0;
291 }
292 /*@
293     MatMultAdd -  Computes v3 = v2 + A * v1.
294 
295   Input Parameters:
296 .    mat - the matrix
297 .    v1, v2 - the vectors
298 
299   Output Parameters:
300 .    v3 - the result
301 
302 .keywords: matrix, multiply, matrix-vector product, add
303 
304 .seealso: MatMultTrans(), MatMult(), MatMultTransAdd()
305 @*/
306 int MatMultAdd(Mat mat,Vec v1,Vec v2,Vec v3)
307 {
308   int ierr;
309   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);PETSCVALIDHEADERSPECIFIC(v1,VEC_COOKIE);
310   PETSCVALIDHEADERSPECIFIC(v2,VEC_COOKIE); PETSCVALIDHEADERSPECIFIC(v3,VEC_COOKIE);
311   PLogEventBegin(MAT_MultAdd,mat,v1,v2,v3);
312   ierr = (*mat->ops.multadd)(mat,v1,v2,v3); CHKERRQ(ierr);
313   PLogEventEnd(MAT_MultAdd,mat,v1,v2,v3);
314   return 0;
315 }
316 /*@
317     MatMultTransAdd - Computes v3 = v2 + A' * v1.
318 
319   Input Parameters:
320 .    mat - the matrix
321 .    v1, v2 - the vectors
322 
323   Output Parameters:
324 .    v3 - the result
325 
326 .keywords: matrix, multiply, matrix-vector product, transpose, add
327 
328 .seealso: MatMultTrans(), MatMultAdd(), MatMult()
329 @*/
330 int MatMultTransAdd(Mat mat,Vec v1,Vec v2,Vec v3)
331 {
332   int ierr;
333   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); PETSCVALIDHEADERSPECIFIC(v1,VEC_COOKIE);
334   PETSCVALIDHEADERSPECIFIC(v2,VEC_COOKIE); PETSCVALIDHEADERSPECIFIC(v3,VEC_COOKIE);
335   if (!mat->ops.multtransadd) SETERRQ(PETSC_ERR_SUP,"MatMultTransAdd");
336   PLogEventBegin(MAT_MultTransAdd,mat,v1,v2,v3);
337   ierr = (*mat->ops.multtransadd)(mat,v1,v2,v3); CHKERRQ(ierr);
338   PLogEventEnd(MAT_MultTransAdd,mat,v1,v2,v3);
339   return 0;
340 }
341 /* ------------------------------------------------------------*/
342 /*@
343    MatGetInfo - Returns information about matrix storage (number of
344    nonzeros, memory).
345 
346    Input Parameters:
347 .  mat - the matrix
348 
349    Output Parameters:
350 .  flag - flag indicating the type of parameters to be returned
351 $    flag = MAT_LOCAL: local matrix
352 $    flag = MAT_GLOBAL_MAX: maximum over all processors
353 $    flag = MAT_GLOBAL_SUM: sum over all processors
354 .   nz - the number of nonzeros
355 .   nzalloc - the number of allocated nonzeros
356 .   mem - the memory used (in bytes)
357 
358 .keywords: matrix, get, info, storage, nonzeros, memory
359 @*/
360 int MatGetInfo(Mat mat,MatInfoType flag,int *nz,int *nzalloc,int *mem)
361 {
362   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
363   if (!mat->ops.getinfo) SETERRQ(PETSC_ERR_SUP,"MatGetInfo");
364   return  (*mat->ops.getinfo)(mat,flag,nz,nzalloc,mem);
365 }
366 /* ----------------------------------------------------------*/
367 /*@
368    MatLUFactor - Performs in-place LU factorization of matrix.
369 
370    Input Parameters:
371 .  mat - the matrix
372 .  row - row permutation
373 .  col - column permutation
374 .  f - expected fill as ratio of original fill.
375 
376 .keywords: matrix, factor, LU, in-place
377 
378 .seealso: MatLUFactorSymbolic(), MatLUFactorNumeric(), MatCholeskyFactor()
379 @*/
380 int MatLUFactor(Mat mat,IS row,IS col,double f)
381 {
382   int ierr;
383   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
384   if (!mat->ops.lufactor) SETERRQ(PETSC_ERR_SUP,"MatLUFactor");
385   PLogEventBegin(MAT_LUFactor,mat,row,col,0);
386   ierr = (*mat->ops.lufactor)(mat,row,col,f); CHKERRQ(ierr);
387   PLogEventEnd(MAT_LUFactor,mat,row,col,0);
388   return 0;
389 }
390 /*@
391    MatLUFactorSymbolic - Performs symbolic LU factorization of matrix.
392    Call this routine before calling MatLUFactorNumeric().
393 
394    Input Parameters:
395 .  mat - the matrix
396 .  row, col - row and column permutations
397 .  f - expected fill as ratio of the original number of nonzeros,
398        for example 3.0; choosing this parameter well can result in
399        more efficient use of time and space.
400 
401    Output Parameters:
402 .  fact - new matrix that has been symbolically factored
403 
404 .keywords: matrix, factor, LU, symbol CHKERRQ(ierr);ic
405 
406 .seealso: MatLUFactor(), MatLUFactorNumeric(), MatCholeskyFactor()
407 @*/
408 int MatLUFactorSymbolic(Mat mat,IS row,IS col,double f,Mat *fact)
409 {
410   int ierr;
411   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
412   if (!fact) SETERRQ(1,"MatLUFactorSymbolic:Missing factor matrix argument");
413   if (!mat->ops.lufactorsymbolic) SETERRQ(PETSC_ERR_SUP,"MatLUFactorSymbolic");
414   OptionsGetDouble(0,"-mat_lu_fill",&f);
415   PLogEventBegin(MAT_LUFactorSymbolic,mat,row,col,0);
416   ierr = (*mat->ops.lufactorsymbolic)(mat,row,col,f,fact); CHKERRQ(ierr);
417   PLogEventEnd(MAT_LUFactorSymbolic,mat,row,col,0);
418   return 0;
419 }
420 /*@
421    MatLUFactorNumeric - Performs numeric LU factorization of a matrix.
422    Call this routine after first calling MatLUFactorSymbolic().
423 
424    Input Parameters:
425 .  mat - the matrix
426 .  row, col - row and  column permutations
427 
428    Output Parameters:
429 .  fact - symbolically factored matrix that must have been generated
430           by MatLUFactorSymbolic()
431 
432    Notes:
433    See MatLUFactor() for in-place factorization.  See
434    MatCholeskyFactorNumeric() for the symmetric, positive definite case.
435 
436 .keywords: matrix, factor, LU, numeric
437 
438 .seealso: MatLUFactorSymbolic(), MatLUFactor(), MatCholeskyFactor()
439 @*/
440 int MatLUFactorNumeric(Mat mat,Mat *fact)
441 {
442   int ierr;
443   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
444   if (!fact) SETERRQ(1,"MatLUFactorNumeric:Missing factor matrix argument");
445   if (!mat->ops.lufactornumeric) SETERRQ(PETSC_ERR_SUP,"MatLUFactorNumeric");
446   PLogEventBegin(MAT_LUFactorNumeric,mat,*fact,0,0);
447   ierr = (*mat->ops.lufactornumeric)(mat,fact); CHKERRQ(ierr);
448   PLogEventEnd(MAT_LUFactorNumeric,mat,*fact,0,0);
449   return 0;
450 }
451 /*@
452    MatCholeskyFactor - Performs in-place Cholesky factorization of a
453    symmetric matrix.
454 
455    Input Parameters:
456 .  mat - the matrix
457 .  perm - row and column permutations
458 .  f - expected fill as ratio of original fill
459 
460    Notes:
461    See MatLUFactor() for the nonsymmetric case.  See also
462    MatCholeskyFactorSymbolic(), and MatCholeskyFactorNumeric().
463 
464 .keywords: matrix, factor, in-place, Cholesky
465 
466 .seealso: MatLUFactor(), MatCholeskyFactorSymbolic()
467 .seealso: MatCholeskyFactorNumeric()
468 @*/
469 int MatCholeskyFactor(Mat mat,IS perm,double f)
470 {
471   int ierr;
472   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
473   if (!mat->ops.choleskyfactor) SETERRQ(PETSC_ERR_SUP,"MatCholeskyFactor");
474   PLogEventBegin(MAT_CholeskyFactor,mat,perm,0,0);
475   ierr = (*mat->ops.choleskyfactor)(mat,perm,f); CHKERRQ(ierr);
476   PLogEventEnd(MAT_CholeskyFactor,mat,perm,0,0);
477   return 0;
478 }
479 /*@
480    MatCholeskyFactorSymbolic - Performs symbolic Cholesky factorization
481    of a symmetric matrix.
482 
483    Input Parameters:
484 .  mat - the matrix
485 .  perm - row and column permutations
486 .  f - expected fill as ratio of original
487 
488    Output Parameter:
489 .  fact - the factored matrix
490 
491    Notes:
492    See MatLUFactorSymbolic() for the nonsymmetric case.  See also
493    MatCholeskyFactor() and MatCholeskyFactorNumeric().
494 
495 .keywords: matrix, factor, factorization, symbolic, Cholesky
496 
497 .seealso: MatLUFactorSymbolic()
498 .seealso: MatCholeskyFactor(), MatCholeskyFactorNumeric()
499 @*/
500 int MatCholeskyFactorSymbolic(Mat mat,IS perm,double f,Mat *fact)
501 {
502   int ierr;
503   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
504   if (!fact)
505     SETERRQ(1,"MatCholeskyFactorSymbolic:Missing factor matrix argument");
506   if (!mat->ops.choleskyfactorsymbolic)SETERRQ(PETSC_ERR_SUP,"MatCholeskyFactorSymbolic");
507   PLogEventBegin(MAT_CholeskyFactorSymbolic,mat,perm,0,0);
508   ierr = (*mat->ops.choleskyfactorsymbolic)(mat,perm,f,fact); CHKERRQ(ierr);
509   PLogEventEnd(MAT_CholeskyFactorSymbolic,mat,perm,0,0);
510   return 0;
511 }
512 /*@
513    MatCholeskyFactorNumeric - Performs numeric Cholesky factorization
514    of a symmetric matrix. Call this routine after first calling
515    MatCholeskyFactorSymbolic().
516 
517    Input Parameter:
518 .  mat - the initial matrix
519 
520    Output Parameter:
521 .  fact - the factored matrix
522 
523 .keywords: matrix, factor, numeric, Cholesky
524 
525 .seealso: MatCholeskyFactorSymbolic(), MatCholeskyFactor()
526 .seealso: MatLUFactorNumeric()
527 @*/
528 int MatCholeskyFactorNumeric(Mat mat,Mat *fact)
529 {
530   int ierr;
531   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
532   if (!fact)
533     SETERRQ(1,"MatCholeskyFactorNumeric:Missing factor matrix argument");
534   if (!mat->ops.choleskyfactornumeric) SETERRQ(PETSC_ERR_SUP,"MatCholeskyFactorNumeric");
535   PLogEventBegin(MAT_CholeskyFactorNumeric,mat,*fact,0,0);
536   ierr = (*mat->ops.choleskyfactornumeric)(mat,fact); CHKERRQ(ierr);
537   PLogEventEnd(MAT_CholeskyFactorNumeric,mat,*fact,0,0);
538   return 0;
539 }
540 /* ----------------------------------------------------------------*/
541 /*@
542    MatSolve - Solves A x = b, given a factored matrix.
543 
544    Input Parameters:
545 .  mat - the factored matrix
546 .  b - the right-hand-side vector
547 
548    Output Parameter:
549 .  x - the result vector
550 
551 .keywords: matrix, linear system, solve, LU, Cholesky, triangular solve
552 
553 .seealso: MatSolveAdd(), MatSolveTrans(), MatSolveTransAdd()
554 @*/
555 int MatSolve(Mat mat,Vec b,Vec x)
556 {
557   int ierr;
558   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
559   PETSCVALIDHEADERSPECIFIC(b,VEC_COOKIE);  PETSCVALIDHEADERSPECIFIC(x,VEC_COOKIE);
560   if (!mat->factor) SETERRQ(1,"MatSolve:Unfactored matrix");
561   if (!mat->ops.solve) SETERRQ(PETSC_ERR_SUP,"MatSolve");
562   PLogEventBegin(MAT_Solve,mat,b,x,0);
563   ierr = (*mat->ops.solve)(mat,b,x); CHKERRQ(ierr);
564   PLogEventEnd(MAT_Solve,mat,b,x,0);
565   return 0;
566 }
567 
568 /* @
569    MatForwardSolve - Solves L x = b, given a factored matrix, A = LU.
570 
571    Input Parameters:
572 .  mat - the factored matrix
573 .  b - the right-hand-side vector
574 
575    Output Parameter:
576 .  x - the result vector
577 
578    Notes:
579    MatSolve() should be used for most applications, as it performs
580    a forward solve followed by a backward solve.
581 
582 .keywords: matrix, forward, LU, Cholesky, triangular solve
583 
584 .seealso: MatSolve(), MatBackwardSolve()
585 @ */
586 int MatForwardSolve(Mat mat,Vec b,Vec x)
587 {
588   int ierr;
589   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
590   PETSCVALIDHEADERSPECIFIC(b,VEC_COOKIE);  PETSCVALIDHEADERSPECIFIC(x,VEC_COOKIE);
591   if (!mat->factor) SETERRQ(1,"MatForwardSolve:Unfactored matrix");
592   if (!mat->ops.forwardsolve) SETERRQ(PETSC_ERR_SUP,"MatForwardSolve");
593   PLogEventBegin(MAT_ForwardSolve,mat,b,x,0);
594   ierr = (*mat->ops.forwardsolve)(mat,b,x); CHKERRQ(ierr);
595   PLogEventEnd(MAT_ForwardSolve,mat,b,x,0);
596   return 0;
597 }
598 
599 /* @
600    MatBackwardSolve - Solves U x = b, given a factored matrix, A = LU.
601 
602    Input Parameters:
603 .  mat - the factored matrix
604 .  b - the right-hand-side vector
605 
606    Output Parameter:
607 .  x - the result vector
608 
609    Notes:
610    MatSolve() should be used for most applications, as it performs
611    a forward solve followed by a backward solve.
612 
613 .keywords: matrix, backward, LU, Cholesky, triangular solve
614 
615 .seealso: MatSolve(), MatForwardSolve()
616 @ */
617 int MatBackwardSolve(Mat mat,Vec b,Vec x)
618 {
619   int ierr;
620   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
621   PETSCVALIDHEADERSPECIFIC(b,VEC_COOKIE);  PETSCVALIDHEADERSPECIFIC(x,VEC_COOKIE);
622   if (!mat->factor) SETERRQ(1,"MatBackwardSolve:Unfactored matrix");
623   if (!mat->ops.backwardsolve) SETERRQ(PETSC_ERR_SUP,"MatBackwardSolve");
624   PLogEventBegin(MAT_BackwardSolve,mat,b,x,0);
625   ierr = (*mat->ops.backwardsolve)(mat,b,x); CHKERRQ(ierr);
626   PLogEventEnd(MAT_BackwardSolve,mat,b,x,0);
627   return 0;
628 }
629 
630 /*@
631    MatSolveAdd - Computes x = y + inv(A)*b, given a factored matrix.
632 
633    Input Parameters:
634 .  mat - the factored matrix
635 .  b - the right-hand-side vector
636 .  y - the vector to be added to
637 
638    Output Parameter:
639 .  x - the result vector
640 
641 .keywords: matrix, linear system, solve, LU, Cholesky, add
642 
643 .seealso: MatSolve(), MatSolveTrans(), MatSolveTransAdd()
644 @*/
645 int MatSolveAdd(Mat mat,Vec b,Vec y,Vec x)
646 {
647   Scalar one = 1.0;
648   Vec    tmp;
649   int    ierr;
650   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);PETSCVALIDHEADERSPECIFIC(y,VEC_COOKIE);
651   PETSCVALIDHEADERSPECIFIC(b,VEC_COOKIE);  PETSCVALIDHEADERSPECIFIC(x,VEC_COOKIE);
652   if (!mat->factor) SETERRQ(1,"MatSolveAdd:Unfactored matrix");
653   PLogEventBegin(MAT_SolveAdd,mat,b,x,y);
654   if (mat->ops.solveadd)  {
655     ierr = (*mat->ops.solveadd)(mat,b,y,x); CHKERRQ(ierr);
656   }
657   else {
658     /* do the solve then the add manually */
659     if (x != y) {
660       ierr = MatSolve(mat,b,x); CHKERRQ(ierr);
661       ierr = VecAXPY(&one,y,x); CHKERRQ(ierr);
662     }
663     else {
664       ierr = VecDuplicate(x,&tmp); CHKERRQ(ierr);
665       PLogObjectParent(mat,tmp);
666       ierr = VecCopy(x,tmp); CHKERRQ(ierr);
667       ierr = MatSolve(mat,b,x); CHKERRQ(ierr);
668       ierr = VecAXPY(&one,tmp,x); CHKERRQ(ierr);
669       ierr = VecDestroy(tmp); CHKERRQ(ierr);
670     }
671   }
672   PLogEventEnd(MAT_SolveAdd,mat,b,x,y);
673   return 0;
674 }
675 /*@
676    MatSolveTrans - Solves A' x = b, given a factored matrix.
677 
678    Input Parameters:
679 .  mat - the factored matrix
680 .  b - the right-hand-side vector
681 
682    Output Parameter:
683 .  x - the result vector
684 
685 .keywords: matrix, linear system, solve, LU, Cholesky, transpose
686 
687 .seealso: MatSolve(), MatSolveAdd(), MatSolveTransAdd()
688 @*/
689 int MatSolveTrans(Mat mat,Vec b,Vec x)
690 {
691   int ierr;
692   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
693   PETSCVALIDHEADERSPECIFIC(b,VEC_COOKIE);  PETSCVALIDHEADERSPECIFIC(x,VEC_COOKIE);
694   if (!mat->factor) SETERRQ(1,"MatSolveTrans:Unfactored matrix");
695   if (!mat->ops.solvetrans) SETERRQ(PETSC_ERR_SUP,"MatSolveTrans");
696   PLogEventBegin(MAT_SolveTrans,mat,b,x,0);
697   ierr = (*mat->ops.solvetrans)(mat,b,x); CHKERRQ(ierr);
698   PLogEventEnd(MAT_SolveTrans,mat,b,x,0);
699   return 0;
700 }
701 /*@
702    MatSolveTransAdd - Computes x = y + inv(trans(A)) b, given a
703                       factored matrix.
704 
705    Input Parameters:
706 .  mat - the factored matrix
707 .  b - the right-hand-side vector
708 .  y - the vector to be added to
709 
710    Output Parameter:
711 .  x - the result vector
712 
713 .keywords: matrix, linear system, solve, LU, Cholesky, transpose, add
714 
715 .seealso: MatSolve(), MatSolveAdd(), MatSolveTrans()
716 @*/
717 int MatSolveTransAdd(Mat mat,Vec b,Vec y,Vec x)
718 {
719   Scalar one = 1.0;
720   int    ierr;
721   Vec    tmp;
722   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);PETSCVALIDHEADERSPECIFIC(y,VEC_COOKIE);
723   PETSCVALIDHEADERSPECIFIC(b,VEC_COOKIE);  PETSCVALIDHEADERSPECIFIC(x,VEC_COOKIE);
724   if (!mat->factor) SETERRQ(1,"MatSolveTransAdd:Unfactored matrix");
725   PLogEventBegin(MAT_SolveTransAdd,mat,b,x,y);
726   if (mat->ops.solvetransadd) {
727     ierr = (*mat->ops.solvetransadd)(mat,b,y,x); CHKERRQ(ierr);
728   }
729   else {
730     /* do the solve then the add manually */
731     if (x != y) {
732       ierr = MatSolveTrans(mat,b,x); CHKERRQ(ierr);
733       ierr = VecAXPY(&one,y,x); CHKERRQ(ierr);
734     }
735     else {
736       ierr = VecDuplicate(x,&tmp); CHKERRQ(ierr);
737       PLogObjectParent(mat,tmp);
738       ierr = VecCopy(x,tmp); CHKERRQ(ierr);
739       ierr = MatSolveTrans(mat,b,x); CHKERRQ(ierr);
740       ierr = VecAXPY(&one,tmp,x); CHKERRQ(ierr);
741       ierr = VecDestroy(tmp); CHKERRQ(ierr);
742     }
743   }
744   PLogEventEnd(MAT_SolveTransAdd,mat,b,x,y);
745   return 0;
746 }
747 /* ----------------------------------------------------------------*/
748 
749 /*@
750    MatRelax - Computes one relaxation sweep.
751 
752    Input Parameters:
753 .  mat - the matrix
754 .  b - the right hand side
755 .  omega - the relaxation factor
756 .  flag - flag indicating the type of SOR, one of
757 $     SOR_FORWARD_SWEEP
758 $     SOR_BACKWARD_SWEEP
759 $     SOR_SYMMETRIC_SWEEP (SSOR method)
760 $     SOR_LOCAL_FORWARD_SWEEP
761 $     SOR_LOCAL_BACKWARD_SWEEP
762 $     SOR_LOCAL_SYMMETRIC_SWEEP (local SSOR)
763 $     SOR_APPLY_UPPER, SOR_APPLY_LOWER - applies
764 $       upper/lower triangular part of matrix to
765 $       vector (with omega)
766 $     SOR_ZERO_INITIAL_GUESS - zero initial guess
767 .  shift -  diagonal shift
768 .  its - the number of iterations
769 
770    Output Parameters:
771 .  x - the solution (can contain an initial guess)
772 
773    Notes:
774    SOR_LOCAL_FORWARD_SWEEP, SOR_LOCAL_BACKWARD_SWEEP, and
775    SOR_LOCAL_SYMMETRIC_SWEEP perform seperate independent smoothings
776    on each processor.
777 
778    Application programmers will not generally use MatRelax() directly,
779    but instead will employ the SLES/PC interface.
780 
781    Notes for Advanced Users:
782    The flags are implemented as bitwise inclusive or operations.
783    For example, use (SOR_ZERO_INITIAL_GUESS | SOR_SYMMETRIC_SWEEP)
784    to specify a zero initial guess for SSOR.
785 
786 .keywords: matrix, relax, relaxation, sweep
787 @*/
788 int MatRelax(Mat mat,Vec b,double omega,MatSORType flag,double shift,
789              int its,Vec x)
790 {
791   int ierr;
792   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
793   PETSCVALIDHEADERSPECIFIC(b,VEC_COOKIE);  PETSCVALIDHEADERSPECIFIC(x,VEC_COOKIE);
794   if (!mat->ops.relax) SETERRQ(PETSC_ERR_SUP,"MatRelax");
795   PLogEventBegin(MAT_Relax,mat,b,x,0);
796   ierr =(*mat->ops.relax)(mat,b,omega,flag,shift,its,x); CHKERRQ(ierr);
797   PLogEventEnd(MAT_Relax,mat,b,x,0);
798   return 0;
799 }
800 
801 /*@C
802    MatConvert - Converts a matrix to another matrix, either of the same
803    or different type.
804 
805    Input Parameters:
806 .  mat - the matrix
807 .  newtype - new matrix type.  Use MATSAME to create a new matrix of the
808    same type as the original matrix.
809 
810    Output Parameter:
811 .  M - pointer to place new matrix
812 
813 .keywords: matrix, copy, convert
814 @*/
815 int MatConvert(Mat mat,MatType newtype,Mat *M)
816 {
817   int ierr;
818   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
819   if (!M) SETERRQ(1,"MatConvert:Bad address");
820   if (newtype == mat->type || newtype == MATSAME) {
821     if (mat->ops.copyprivate) {
822       PLogEventBegin(MAT_Convert,mat,0,0,0);
823       ierr = (*mat->ops.copyprivate)(mat,M); CHKERRQ(ierr);
824       PLogEventEnd(MAT_Convert,mat,0,0,0);
825       return 0;
826     }
827   }
828   if (!mat->ops.convert) SETERRQ(PETSC_ERR_SUP,"MatConvert");
829   PLogEventBegin(MAT_Convert,mat,0,0,0);
830   if (mat->ops.convert) {
831     ierr = (*mat->ops.convert)(mat,newtype,M); CHKERRQ(ierr);
832   }
833   PLogEventEnd(MAT_Convert,mat,0,0,0);
834   return 0;
835 }
836 
837 /*@
838    MatGetDiagonal - Gets the diagonal of a matrix.
839 
840    Input Parameters:
841 .  mat - the matrix
842 
843    Output Parameters:
844 .  v - the vector for storing the diagonal
845 
846 .keywords: matrix, get, diagonal
847 @*/
848 int MatGetDiagonal(Mat mat,Vec v)
849 {
850   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
851   PETSCVALIDHEADERSPECIFIC(v,VEC_COOKIE);
852   if (mat->ops.getdiagonal) return (*mat->ops.getdiagonal)(mat,v);
853   SETERRQ(PETSC_ERR_SUP,"MatGetDiagonal");
854 }
855 
856 /*@
857    MatTranspose - Computes an in-place or out-of-place transpose of a matrix.
858 
859    Input Parameters:
860 .  mat - the matrix to transpose
861 
862    Output Parameters:
863 .   B - the transpose -  pass in zero for an in-place transpose
864 
865 .keywords: matrix, transpose
866 @*/
867 int MatTranspose(Mat mat,Mat *B)
868 {
869   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
870   if (mat->ops.transpose) return (*mat->ops.transpose)(mat,B);
871   SETERRQ(PETSC_ERR_SUP,"MatTranspose");
872 }
873 
874 /*@
875    MatEqual - Compares two matrices.  Returns 1 if two matrices are equal.
876 
877    Input Parameters:
878 .  mat1 - the first matrix
879 .  mat2 - the second matrix
880 
881    Returns:
882    Returns 1 if the matrices are equal; returns 0 otherwise.
883 
884 .keywords: matrix, equal, equivalent
885 @*/
886 int MatEqual(Mat mat1,Mat mat2)
887 {
888   PETSCVALIDHEADERSPECIFIC(mat1,MAT_COOKIE); PETSCVALIDHEADERSPECIFIC(mat2,MAT_COOKIE);
889   if (mat1->ops.equal) return (*mat1->ops.equal)(mat1,mat2);
890   SETERRQ(PETSC_ERR_SUP,"MatEqual");
891 }
892 
893 /*@
894    MatScale - Scales a matrix on the left and right by diagonal
895    matrices that are stored as vectors.  Either of the two scaling
896    matrices can be null.
897 
898    Input Parameters:
899 .  mat - the matrix to be scaled
900 .  l - the left scaling vector
901 .  r - the right scaling vector
902 
903 .keywords: matrix, scale
904 @*/
905 int MatScale(Mat mat,Vec l,Vec r)
906 {
907   int ierr;
908   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
909   if (!mat->ops.scale) SETERRQ(PETSC_ERR_SUP,"MatScale");
910   if (l) PETSCVALIDHEADERSPECIFIC(l,VEC_COOKIE);
911   if (r) PETSCVALIDHEADERSPECIFIC(r,VEC_COOKIE);
912   PLogEventBegin(MAT_Scale,mat,0,0,0);
913   ierr = (*mat->ops.scale)(mat,l,r); CHKERRQ(ierr);
914   PLogEventEnd(MAT_Scale,mat,0,0,0);
915   return 0;
916 }
917 
918 /*@
919    MatNorm - Calculates various norms of a matrix.
920 
921    Input Parameters:
922 .  mat - the matrix
923 .  type - the type of norm, NORM_1, NORM_2, NORM_FROBENIUS, NORM_INFINITY
924 
925    Output Parameters:
926 .  norm - the resulting norm
927 
928 .keywords: matrix, norm, Frobenius
929 @*/
930 int MatNorm(Mat mat,MatNormType type,double *norm)
931 {
932   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
933   if (!norm) SETERRQ(1,"MatNorm:bad addess for value");
934   if (mat->ops.norm) return (*mat->ops.norm)(mat,type,norm);
935   SETERRQ(PETSC_ERR_SUP,"MatNorm:Not for this matrix type");
936 }
937 
938 /*@
939    MatAssemblyBegin - Begins assembling the matrix.  This routine should
940    be called after completing all calls to MatSetValues().
941 
942    Input Parameters:
943 .  mat - the matrix
944 .  type - type of assembly, either FLUSH_ASSEMBLY or FINAL_ASSEMBLY
945 
946    Notes:
947    MatSetValues() generally caches the values.  The matrix is ready to
948    use only after MatAssemblyBegin() and MatAssemblyEnd() have been called.
949    Use FLUSH_ASSEMBLY when switching between ADD_VALUES and SetValues; use
950    FINAL_ASSEMBLY for the final assembly before the matrix is used.
951 
952 .keywords: matrix, assembly, assemble, begin
953 
954 .seealso: MatAssemblyEnd(), MatSetValues()
955 @*/
956 int MatAssemblyBegin(Mat mat,MatAssemblyType type)
957 {
958   int ierr;
959   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
960   PLogEventBegin(MAT_AssemblyBegin,mat,0,0,0);
961   if (mat->ops.assemblybegin) {ierr = (*mat->ops.assemblybegin)(mat,type); CHKERRQ(ierr);}
962   PLogEventEnd(MAT_AssemblyBegin,mat,0,0,0);
963   return 0;
964 }
965 /*@
966    MatAssemblyEnd - Completes assembling the matrix.  This routine should
967    be called after all calls to MatSetValues() and after MatAssemblyBegin().
968 
969    Input Parameters:
970 .  mat - the matrix
971 .  type - type of assembly, either FLUSH_ASSEMBLY or FINAL_ASSEMBLY
972 
973    Note:
974    MatSetValues() generally caches the values.  The matrix is ready to
975    use only after MatAssemblyBegin() and MatAssemblyEnd() have been called.
976    Use FLUSH_ASSEMBLY when switching between ADD_VALUES and SetValues; use
977    FINAL_ASSEMBLY for the final assembly before the matrix is used.
978 
979 .keywords: matrix, assembly, assemble, end
980 
981 .seealso: MatAssemblyBegin(), MatSetValues()
982 @*/
983 int MatAssemblyEnd(Mat mat,MatAssemblyType type)
984 {
985   int ierr;
986   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
987   PLogEventBegin(MAT_AssemblyEnd,mat,0,0,0);
988   if (mat->ops.assemblyend) {ierr = (*mat->ops.assemblyend)(mat,type); CHKERRQ(ierr);}
989   PLogEventEnd(MAT_AssemblyEnd,mat,0,0,0);
990   return 0;
991 }
992 /*@
993    MatCompress - Tries to store the matrix in as little space as
994    possible.  May fail if memory is already fully used, since it
995    tries to allocate new space.
996 
997    Input Parameters:
998 .  mat - the matrix
999 
1000 .keywords: matrix, compress
1001 @*/
1002 int MatCompress(Mat mat)
1003 {
1004   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
1005   if (mat->ops.compress) return (*mat->ops.compress)(mat);
1006   return 0;
1007 }
1008 /*@
1009    MatSetOption - Sets a parameter option for a matrix. Some options
1010    may be specific to certain storage formats.  Some options
1011    determine how values will be inserted (or added). Sorted,
1012    row-oriented input will generally assemble the fastest. The default
1013    is row-oriented, nonsorted input.
1014 
1015    Input Parameters:
1016 .  mat - the matrix
1017 .  option - the option, one of the following:
1018 $    ROW_ORIENTED, COLUMN_ORIENTED, ROWS_SORTED,
1019 $    COLUMNS_SORTED, NO_NEW_NONZERO_LOCATIONS,
1020 $    SYMMETRIC_MATRIX,
1021 $    YES_NEW_NONZERO_LOCATIONS, and possibly others
1022 
1023    Notes:
1024    If using a Fortran 77 module to compute a matrix, one may need to
1025    use the column-oriented option (or convert to the row-oriented
1026    format).
1027 
1028    NO_NEW_NONZERO_LOCATIONS indicates that any add or insertion
1029    that will generate a new entry in the nonzero structure is ignored.
1030    What this means is if memory is not allocated for this particular
1031    lot, then the insertion is ignored. For dense matrices where
1032    the entire array is allocated no entries are ever ignored. This
1033    may not be a good idea??
1034 
1035 .keywords: matrix, option, row-oriented, column-oriented, sorted, nonzero
1036 @*/
1037 int MatSetOption(Mat mat,MatOption op)
1038 {
1039   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
1040   if (mat->ops.setoption) return (*mat->ops.setoption)(mat,op);
1041   return 0;
1042 }
1043 
1044 /*@
1045    MatZeroEntries - Zeros all entries of a matrix.  For sparse matrices
1046    this routine retains the old nonzero structure.
1047 
1048    Input Parameters:
1049 .  mat - the matrix
1050 
1051 .keywords: matrix, zero, entries
1052 
1053 .seealso: MatZeroRows()
1054 @*/
1055 int MatZeroEntries(Mat mat)
1056 {
1057   int ierr;
1058   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
1059   if (!mat->ops.zeroentries) SETERRQ(PETSC_ERR_SUP,"MatZeroEntries");
1060   PLogEventBegin(MAT_ZeroEntries,mat,0,0,0);
1061   ierr = (*mat->ops.zeroentries)(mat); CHKERRQ(ierr);
1062   PLogEventEnd(MAT_ZeroEntries,mat,0,0,0);
1063   return 0;
1064 }
1065 
1066 /*@
1067    MatZeroRows - Zeros all entries (except possibly the main diagonal)
1068    of a set of rows of a matrix.
1069 
1070    Input Parameters:
1071 .  mat - the matrix
1072 .  is - index set of rows to remove
1073 .  diag - pointer to value put in all diagonals of eliminated rows.
1074           Note that diag is not a pointer to an array, but merely a
1075           pointer to a single value.
1076 
1077    Notes:
1078    For the AIJ and row matrix formats this removes the old nonzero
1079    structure, but does not release memory.  For the dense and block
1080    diagonal formats this does not alter the nonzero structure.
1081 
1082    The user can set a value in the diagonal entry (or for the AIJ and
1083    row formats can optionally remove the main diagonal entry from the
1084    nonzero structure as well, by passing a null pointer as the final
1085    argument).
1086 
1087 .keywords: matrix, zero, rows, boundary conditions
1088 
1089 .seealso: MatZeroEntries(), MatGetSubMatrix(), MatGetSubMatrixInPlace()
1090 @*/
1091 int MatZeroRows(Mat mat,IS is, Scalar *diag)
1092 {
1093   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
1094   if (mat->ops.zerorows) return (*mat->ops.zerorows)(mat,is,diag);
1095   SETERRQ(PETSC_ERR_SUP,"MatZeroRows");
1096 }
1097 
1098 /*@
1099    MatGetSize - Returns the numbers of rows and columns in a matrix.
1100 
1101    Input Parameter:
1102 .  mat - the matrix
1103 
1104    Output Parameters:
1105 .  m - the number of global rows
1106 .  n - the number of global columns
1107 
1108 .keywords: matrix, dimension, size, rows, columns, global, get
1109 
1110 .seealso: MatGetLocalSize()
1111 @*/
1112 int MatGetSize(Mat mat,int *m,int* n)
1113 {
1114   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
1115   if (!m || !n) SETERRQ(1,"MatGetSize:Bad address for result");
1116   return (*mat->ops.getsize)(mat,m,n);
1117 }
1118 
1119 /*@
1120    MatGetLocalSize - Returns the number of rows and columns in a matrix
1121    stored locally.  This information may be implementation dependent, so
1122    use with care.
1123 
1124    Input Parameters:
1125 .  mat - the matrix
1126 
1127    Output Parameters:
1128 .  m - the number of local rows
1129 .  n - the number of local columns
1130 
1131 .keywords: matrix, dimension, size, local, rows, columns, get
1132 
1133 .seealso: MatGetSize()
1134 @*/
1135 int MatGetLocalSize(Mat mat,int *m,int* n)
1136 {
1137   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
1138   if (!m || !n) SETERRQ(1,"MatGetLocalSize:Bad address for result");
1139   return (*mat->ops.getlocalsize)(mat,m,n);
1140 }
1141 
1142 /*@
1143    MatGetOwnershipRange - Returns the range of matrix rows owned by
1144    this processor, assuming that the matrix is laid out with the first
1145    n1 rows on the first processor, the next n2 rows on the second, etc.
1146    For certain parallel layouts this range may not be well-defined.
1147 
1148    Input Parameters:
1149 .  mat - the matrix
1150 
1151    Output Parameters:
1152 .  m - the first local row
1153 .  n - one more then the last local row
1154 
1155 .keywords: matrix, get, range, ownership
1156 @*/
1157 int MatGetOwnershipRange(Mat mat,int *m,int* n)
1158 {
1159   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
1160   if (!m || !n) SETERRQ(1,"MatGetOwnershipRange:Bad address for result");
1161   if (mat->ops.getownershiprange) return (*mat->ops.getownershiprange)(mat,m,n);
1162   SETERRQ(PETSC_ERR_SUP,"MatGetOwnershipRange");
1163 }
1164 
1165 /*@
1166    MatILUFactorSymbolic - Performs symbolic ILU factorization of a matrix.
1167    Uses levels of fill only, not drop tolerance. Use MatLUFactorNumeric()
1168    to complete the factorization.
1169 
1170    Input Parameters:
1171 .  mat - the matrix
1172 .  row - row permutation
1173 .  column - column permutation
1174 .  fill - number of levels of fill
1175 .  f - expected fill as ratio of original fill
1176 
1177    Output Parameters:
1178 .  fact - puts factor
1179 
1180 .keywords: matrix, factor, incomplete, ILU, symbolic, fill
1181 
1182 .seealso: MatLUFactorSymbolic(), MatLUFactorNumeric()
1183 @*/
1184 int MatILUFactorSymbolic(Mat mat,IS row,IS col,double f,int fill,Mat *fact)
1185 {
1186   int ierr;
1187   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
1188   if (fill < 0) SETERRQ(1,"MatILUFactorSymbolic:Levels of fill negative");
1189   if (!fact) SETERRQ(1,"MatILUFactorSymbolic:Fact argument is missing");
1190   if (!mat->ops.ilufactorsymbolic) SETERRQ(PETSC_ERR_SUP,"MatILUFactorSymbolic");
1191   PLogEventBegin(MAT_ILUFactorSymbolic,mat,row,col,0);
1192   ierr = (*mat->ops.ilufactorsymbolic)(mat,row,col,f,fill,fact);
1193   CHKERRQ(ierr);
1194   PLogEventEnd(MAT_ILUFactorSymbolic,mat,row,col,0);
1195   return 0;
1196 }
1197 
1198 /*@
1199    MatIncompleteCholeskyFactorSymbolic - Performs symbolic incomplete
1200    Cholesky factorization for a symmetric matrix.  Use
1201    MatCholeskyFactorNumeric() to complete the factorization.
1202 
1203    Input Parameters:
1204 .  mat - the matrix
1205 .  perm - row and column permutation
1206 .  fill - levels of fill
1207 .  f - expected fill as ratio of original fill
1208 
1209    Output Parameter:
1210 .  fact - the factored matrix
1211 
1212    Note:  Currently only no-fill factorization is supported.
1213 
1214 .keywords: matrix, factor, incomplete, ICC, Cholesky, symbolic, fill
1215 
1216 .seealso: MatCholeskyFactorNumeric(), MatCholeskyFactor()
1217 @*/
1218 int MatIncompleteCholeskyFactorSymbolic(Mat mat,IS perm,double f,int fill,
1219                                         Mat *fact)
1220 {
1221   int ierr;
1222   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
1223   if (fill < 0) SETERRQ(1,"MatIncompleteCholeskyFactorSymbolic:Fill negative");
1224   if (!fact) SETERRQ(1,"MatIncompleteCholeskyFactorSymbolic:Missing fact argument");
1225   if (!mat->ops.incompletecholeskyfactorsymbolic)
1226     SETERRQ(PETSC_ERR_SUP,"MatIncompleteCholeskyFactorSymbolic");
1227   PLogEventBegin(MAT_IncompleteCholeskyFactorSymbolic,mat,perm,0,0);
1228   ierr = (*mat->ops.incompletecholeskyfactorsymbolic)(mat,perm,f,fill,fact);
1229   CHKERRQ(ierr);
1230   PLogEventEnd(MAT_IncompleteCholeskyFactorSymbolic,mat,perm,0,0);
1231   return 0;
1232 }
1233 
1234 /*@C
1235    MatGetArray - Returns a pointer to the element values in the matrix.
1236    This routine  is implementation dependent, and may not even work for
1237    certain matrix types.
1238 
1239    Input Parameter:
1240 .  mat - the matrix
1241 
1242    Output Parameter:
1243 .  v - the location of the values
1244 
1245 .keywords: matrix, array, elements, values
1246 @*/
1247 int MatGetArray(Mat mat,Scalar **v)
1248 {
1249   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
1250   if (!v) SETERRQ(1,"MatGetArray:Bad input, array pointer location");
1251   if (!mat->ops.getarray) SETERRQ(PETSC_ERR_SUP,"MatGetArraye");
1252   return (*mat->ops.getarray)(mat,v);
1253 }
1254 
1255 /*@
1256    MatGetSubMatrix - Extracts a submatrix from a matrix.
1257 
1258    Input Parameters:
1259 .  mat - the matrix
1260 .  irow, icol - index sets of rows and columns to extract
1261 
1262    Output Parameter:
1263 .  submat - the submatrix
1264 
1265    Notes:
1266    MatGetSubMatrix() can be useful in setting boundary conditions.
1267 
1268 .keywords: matrix, get, submatrix, boundary conditions
1269 
1270 .seealso: MatZeroRows(), MatGetSubMatrixInPlace()
1271 @*/
1272 int MatGetSubMatrix(Mat mat,IS irow,IS icol,Mat *submat)
1273 {
1274   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
1275   if (!mat->ops.getsubmatrix) SETERRQ(PETSC_ERR_SUP,"MatGetSubMatrix");
1276   return (*mat->ops.getsubmatrix)(mat,irow,icol,submat);
1277 }
1278 
1279 /*@
1280    MatGetSubMatrixInPlace - Extracts a submatrix from a matrix, returning
1281    the submatrix in place of the original matrix.
1282 
1283    Input Parameters:
1284 .  mat - the matrix
1285 .  irow, icol - index sets of rows and columns to extract
1286 
1287 .keywords: matrix, get, submatrix, boundary conditions, in-place
1288 
1289 .seealso: MatZeroRows(), MatGetSubMatrix()
1290 @*/
1291 int MatGetSubMatrixInPlace(Mat mat,IS irow,IS icol)
1292 {
1293   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
1294   if (!mat->ops.getsubmatrixinplace) SETERRQ(PETSC_ERR_SUP,"MatGetSubmatrixInPlace");
1295   return (*mat->ops.getsubmatrixinplace)(mat,irow,icol);
1296 }
1297 
1298 /*@
1299    MatGetType - Returns the type of the matrix, one of MATSEQDENSE, MATSEQAIJ, etc.
1300 
1301    Input Parameter:
1302 .  mat - the matrix
1303 
1304    Ouput Parameter:
1305 .  type - the matrix type
1306 
1307    Notes:
1308    The matrix types are given in petsc/include/mat.h
1309 
1310 .keywords: matrix, get, type
1311 
1312 .seealso:  MatGetName()
1313 @*/
1314 int MatGetType(Mat mat,MatType *type)
1315 {
1316   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
1317   *type = (MatType) mat->type;
1318   return 0;
1319 }
1320