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