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