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