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