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