xref: /petsc/src/mat/interface/matrix.c (revision 567e79a179090e2aa3c1c687070390f6bb7fc8f3)
1 #ifndef lint
2 static char vcid[] = "$Id: matrix.c,v 1.106 1995/10/26 19:17:24 curfman Exp curfman $";
3 #endif
4 
5 /*
6    This is where the abstract matrix operations are defined
7 */
8 
9 #include "petsc.h"
10 #include "matimpl.h"        /*I "mat.h" I*/
11 #include "vec/vecimpl.h"
12 #include "pinclude/pviewer.h"
13 #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 MATSEQAIJ:
873         ierr = MatCreateSeqAIJ(mat->comm,m,n,0,0,M); CHKERRQ(ierr);
874         break;
875       case MATSEQROW:
876         ierr = MatCreateSeqRow(mat->comm,m,n,0,0,M); CHKERRQ(ierr);
877         break;
878       case MATMPIROW:
879         ierr = MatCreateMPIRow(MPI_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,
880                                m,n,0,0,0,0,M); CHKERRQ(ierr);
881         break;
882       case MATMPIROWBS:
883         if (m != n) SETERRQ(1,"MatConvert:MATMPIROWBS matrix must be square");
884         ierr = MatCreateMPIRowbs(MPI_COMM_WORLD,PETSC_DECIDE,m,0,0,0,M);
885                CHKERRQ(ierr);
886         break;
887       case MATMPIAIJ:
888         ierr = MatCreateMPIAIJ(MPI_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,
889                                m,n,0,0,0,0,M); CHKERRQ(ierr);
890         break;
891       case MATSEQDENSE:
892         ierr = MatCreateSeqDense(mat->comm,m,n,M); CHKERRQ(ierr);
893         break;
894       case MATMPIDENSE:
895         ierr = MatCreateMPIDense(MPI_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,
896                                  m,n,M); CHKERRQ(ierr);
897         break;
898       case MATSEQBDIAG:
899         {
900         int nb = 1; /* Default block size = 1 */
901         OptionsGetInt(0,"-mat_bdiag_bsize",&nb);
902         ierr = MatCreateSeqBDiag(mat->comm,m,n,0,nb,0,0,M); CHKERRQ(ierr);
903         break;
904         }
905       case MATMPIBDIAG:
906         {
907         int nb = 1; /* Default block size = 1 */
908         OptionsGetInt(0,"-mat_bdiag_bsize",&nb);
909         ierr = MatCreateMPIBDiag(MPI_COMM_WORLD,PETSC_DECIDE,m,n,0,nb,0,0,M);
910         CHKERRQ(ierr);
911         break;
912         }
913       default:
914         SETERRQ(1,"MatConvert:Matrix type is not currently supported");
915     }
916     ierr = MatGetOwnershipRange(*M,&rstart,&rend); CHKERRQ(ierr);
917     for (i=rstart; i<rend; i++) {
918       ierr = MatGetRow(mat,i,&nz,&cwork,&vwork); CHKERRQ(ierr);
919       ierr = MatSetValues(*M,1,&i,nz,cwork,vwork,INSERT_VALUES); CHKERRQ(ierr);
920       ierr = MatRestoreRow(mat,i,&nz,&cwork,&vwork); CHKERRQ(ierr);
921     }
922     ierr = MatAssemblyBegin(*M,FINAL_ASSEMBLY); CHKERRQ(ierr);
923     ierr = MatAssemblyEnd(*M,FINAL_ASSEMBLY); CHKERRQ(ierr);
924   }
925   else {
926     ierr = (*mat->ops.convert)(mat,newtype,M); CHKERRQ(ierr);
927   }
928   PLogEventEnd(MAT_Convert,mat,0,0,0);
929   return 0;
930 }
931 
932 /*@
933    MatGetDiagonal - Gets the diagonal of a matrix.
934 
935    Input Parameters:
936 .  mat - the matrix
937 
938    Output Parameters:
939 .  v - the vector for storing the diagonal
940 
941 .keywords: matrix, get, diagonal
942 @*/
943 int MatGetDiagonal(Mat mat,Vec v)
944 {
945   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
946   PETSCVALIDHEADERSPECIFIC(v,VEC_COOKIE);
947   if (mat->ops.getdiagonal) return (*mat->ops.getdiagonal)(mat,v);
948   SETERRQ(PETSC_ERR_SUP,"MatGetDiagonal");
949 }
950 
951 /*@C
952    MatTranspose - Computes an in-place or out-of-place transpose of a matrix.
953 
954    Input Parameters:
955 .  mat - the matrix to transpose
956 
957    Output Parameters:
958 .   B - the transpose -  pass in zero for an in-place transpose
959 
960 .keywords: matrix, transpose
961 @*/
962 int MatTranspose(Mat mat,Mat *B)
963 {
964   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
965   if (mat->ops.transpose) return (*mat->ops.transpose)(mat,B);
966   SETERRQ(PETSC_ERR_SUP,"MatTranspose");
967 }
968 
969 /*@
970    MatEqual - Compares two matrices.  Returns 1 if two matrices are equal.
971 
972    Input Parameters:
973 .  mat1 - the first matrix
974 .  mat2 - the second matrix
975 
976    Returns:
977    Returns 1 if the matrices are equal; returns 0 otherwise.
978 
979 .keywords: matrix, equal, equivalent
980 @*/
981 int MatEqual(Mat mat1,Mat mat2)
982 {
983   PETSCVALIDHEADERSPECIFIC(mat1,MAT_COOKIE); PETSCVALIDHEADERSPECIFIC(mat2,MAT_COOKIE);
984   if (mat1->ops.equal) return (*mat1->ops.equal)(mat1,mat2);
985   SETERRQ(PETSC_ERR_SUP,"MatEqual");
986 }
987 
988 /*@
989    MatScale - Scales a matrix on the left and right by diagonal
990    matrices that are stored as vectors.  Either of the two scaling
991    matrices can be null.
992 
993    Input Parameters:
994 .  mat - the matrix to be scaled
995 .  l - the left scaling vector
996 .  r - the right scaling vector
997 
998 .keywords: matrix, scale
999 @*/
1000 int MatScale(Mat mat,Vec l,Vec r)
1001 {
1002   int ierr;
1003   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
1004   if (!mat->ops.scale) SETERRQ(PETSC_ERR_SUP,"MatScale");
1005   if (l) PETSCVALIDHEADERSPECIFIC(l,VEC_COOKIE);
1006   if (r) PETSCVALIDHEADERSPECIFIC(r,VEC_COOKIE);
1007   PLogEventBegin(MAT_Scale,mat,0,0,0);
1008   ierr = (*mat->ops.scale)(mat,l,r); CHKERRQ(ierr);
1009   PLogEventEnd(MAT_Scale,mat,0,0,0);
1010   return 0;
1011 }
1012 
1013 /*@
1014    MatNorm - Calculates various norms of a matrix.
1015 
1016    Input Parameters:
1017 .  mat - the matrix
1018 .  type - the type of norm, NORM_1, NORM_2, NORM_FROBENIUS, NORM_INFINITY
1019 
1020    Output Parameters:
1021 .  norm - the resulting norm
1022 
1023 .keywords: matrix, norm, Frobenius
1024 @*/
1025 int MatNorm(Mat mat,MatNormType type,double *norm)
1026 {
1027   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
1028   if (!norm) SETERRQ(1,"MatNorm:bad addess for value");
1029   if (mat->ops.norm) return (*mat->ops.norm)(mat,type,norm);
1030   SETERRQ(PETSC_ERR_SUP,"MatNorm:Not for this matrix type");
1031 }
1032 
1033 /*@
1034    MatAssemblyBegin - Begins assembling the matrix.  This routine should
1035    be called after completing all calls to MatSetValues().
1036 
1037    Input Parameters:
1038 .  mat - the matrix
1039 .  type - type of assembly, either FLUSH_ASSEMBLY or FINAL_ASSEMBLY
1040 
1041    Notes:
1042    MatSetValues() generally caches the values.  The matrix is ready to
1043    use only after MatAssemblyBegin() and MatAssemblyEnd() have been called.
1044    Use FLUSH_ASSEMBLY when switching between ADD_VALUES and SetValues; use
1045    FINAL_ASSEMBLY for the final assembly before the matrix is used.
1046 
1047 .keywords: matrix, assembly, assemble, begin
1048 
1049 .seealso: MatAssemblyEnd(), MatSetValues()
1050 @*/
1051 int MatAssemblyBegin(Mat mat,MatAssemblyType type)
1052 {
1053   int ierr;
1054   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
1055   PLogEventBegin(MAT_AssemblyBegin,mat,0,0,0);
1056   if (mat->ops.assemblybegin) {ierr = (*mat->ops.assemblybegin)(mat,type); CHKERRQ(ierr);}
1057   PLogEventEnd(MAT_AssemblyBegin,mat,0,0,0);
1058   return 0;
1059 }
1060 
1061 /*@
1062    MatAssemblyEnd - Completes assembling the matrix.  This routine should
1063    be called after all calls to MatSetValues() and after MatAssemblyBegin().
1064 
1065    Input Parameters:
1066 .  mat - the matrix
1067 .  type - type of assembly, either FLUSH_ASSEMBLY or FINAL_ASSEMBLY
1068 
1069    Options Database Keys:
1070 $  -mat_view_draw : Draw nonzero structure of matrix at conclusion of MatEndAssembly(),
1071                using MatView() and DrawOpenX().
1072 $  -mat_view_info : Prints info on matrix.
1073 $  -mat_view_info_detailed: More detailed information.
1074 $  -mat_view_ascii : Prints matrix out in ascii.
1075 $  -display <name> : Set display name (default is host)
1076 $  -pause <sec> : Set number of seconds to pause after display
1077 
1078    Note:
1079    MatSetValues() generally caches the values.  The matrix is ready to
1080    use only after MatAssemblyBegin() and MatAssemblyEnd() have been called.
1081    Use FLUSH_ASSEMBLY when switching between ADD_VALUES and SetValues; use
1082    FINAL_ASSEMBLY for the final assembly before the matrix is used.
1083 
1084 .keywords: matrix, assembly, assemble, end
1085 
1086 .seealso: MatAssemblyBegin(), MatSetValues()
1087 @*/
1088 int MatAssemblyEnd(Mat mat,MatAssemblyType type)
1089 {
1090   int        ierr;
1091   static int inassm = 0;
1092   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
1093   inassm++;
1094   PLogEventBegin(MAT_AssemblyEnd,mat,0,0,0);
1095   if (mat->ops.assemblyend) {ierr = (*mat->ops.assemblyend)(mat,type); CHKERRQ(ierr);}
1096   PLogEventEnd(MAT_AssemblyEnd,mat,0,0,0);
1097   if (inassm == 1) {
1098     if (OptionsHasName(0,"-mat_view_info")) {
1099       Viewer viewer;
1100       ierr = ViewerFileOpenASCII(mat->comm,"stdout",&viewer);CHKERRQ(ierr);
1101       ierr = ViewerFileSetFormat(viewer,FILE_FORMAT_INFO,0);CHKERRQ(ierr);
1102       ierr = MatView(mat,viewer); CHKERRQ(ierr);
1103       ierr = ViewerDestroy(viewer); CHKERRQ(ierr);
1104     }
1105     if (OptionsHasName(0,"-mat_view_info_detailed")) {
1106       Viewer viewer;
1107       ierr = ViewerFileOpenASCII(mat->comm,"stdout",&viewer);CHKERRQ(ierr);
1108       ierr = ViewerFileSetFormat(viewer,FILE_FORMAT_INFO_DETAILED,0);CHKERRQ(ierr);
1109       ierr = MatView(mat,viewer); CHKERRQ(ierr);
1110       ierr = ViewerDestroy(viewer); CHKERRQ(ierr);
1111     }
1112     if (OptionsHasName(0,"-mat_view_draw")) {
1113       DrawCtx    win;
1114       ierr = DrawOpenX(mat->comm,0,0,0,0,300,300,&win); CHKERRQ(ierr);
1115       ierr = MatView(mat,(Viewer)win); CHKERRQ(ierr);
1116       ierr = DrawSyncFlush(win); CHKERRQ(ierr);
1117       ierr = DrawDestroy(win); CHKERRQ(ierr);
1118     }
1119   }
1120   inassm--;
1121   return 0;
1122 }
1123 
1124 /*@
1125    MatCompress - Tries to store the matrix in as little space as
1126    possible.  May fail if memory is already fully used, since it
1127    tries to allocate new space.
1128 
1129    Input Parameters:
1130 .  mat - the matrix
1131 
1132 .keywords: matrix, compress
1133 @*/
1134 int MatCompress(Mat mat)
1135 {
1136   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
1137   if (mat->ops.compress) return (*mat->ops.compress)(mat);
1138   return 0;
1139 }
1140 /*@
1141    MatSetOption - Sets a parameter option for a matrix. Some options
1142    may be specific to certain storage formats.  Some options
1143    determine how values will be inserted (or added). Sorted,
1144    row-oriented input will generally assemble the fastest. The default
1145    is row-oriented, nonsorted input.
1146 
1147    Input Parameters:
1148 .  mat - the matrix
1149 .  option - the option, one of the following:
1150 $    ROW_ORIENTED
1151 $    COLUMN_ORIENTED,
1152 $    ROWS_SORTED,
1153 $    COLUMNS_SORTED,
1154 $    NO_NEW_NONZERO_LOCATIONS,
1155 $    YES_NEW_NONZERO_LOCATIONS,
1156 $    SYMMETRIC_MATRIX,
1157 $    STRUCTURALLY_SYMMETRIC_MATRIX,
1158 $    NO_NEW_DIAGONALS,
1159 $    YES_NEW_DIAGONALS,
1160 $    and possibly others.
1161 
1162    Notes:
1163    Some options are relevant only for particular matrix types and
1164    are thus ignored by others.  Other options are not supported by
1165    certain matrix types and will generate an error message if set.
1166 
1167    If using a Fortran 77 module to compute a matrix, one may need to
1168    use the column-oriented option (or convert to the row-oriented
1169    format).
1170 
1171    NO_NEW_NONZERO_LOCATIONS indicates that any add or insertion
1172    that will generate a new entry in the nonzero structure is ignored.
1173    What this means is if memory is not allocated for this particular
1174    lot, then the insertion is ignored. For dense matrices, where
1175    the entire array is allocated, no entries are ever ignored.
1176 
1177 .keywords: matrix, option, row-oriented, column-oriented, sorted, nonzero
1178 @*/
1179 int MatSetOption(Mat mat,MatOption op)
1180 {
1181   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
1182   if (mat->ops.setoption) return (*mat->ops.setoption)(mat,op);
1183   return 0;
1184 }
1185 
1186 /*@
1187    MatZeroEntries - Zeros all entries of a matrix.  For sparse matrices
1188    this routine retains the old nonzero structure.
1189 
1190    Input Parameters:
1191 .  mat - the matrix
1192 
1193 .keywords: matrix, zero, entries
1194 
1195 .seealso: MatZeroRows()
1196 @*/
1197 int MatZeroEntries(Mat mat)
1198 {
1199   int ierr;
1200   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
1201   if (!mat->ops.zeroentries) SETERRQ(PETSC_ERR_SUP,"MatZeroEntries");
1202   PLogEventBegin(MAT_ZeroEntries,mat,0,0,0);
1203   ierr = (*mat->ops.zeroentries)(mat); CHKERRQ(ierr);
1204   PLogEventEnd(MAT_ZeroEntries,mat,0,0,0);
1205   return 0;
1206 }
1207 
1208 /*@
1209    MatZeroRows - Zeros all entries (except possibly the main diagonal)
1210    of a set of rows of a matrix.
1211 
1212    Input Parameters:
1213 .  mat - the matrix
1214 .  is - index set of rows to remove
1215 .  diag - pointer to value put in all diagonals of eliminated rows.
1216           Note that diag is not a pointer to an array, but merely a
1217           pointer to a single value.
1218 
1219    Notes:
1220    For the AIJ and row matrix formats this removes the old nonzero
1221    structure, but does not release memory.  For the dense and block
1222    diagonal formats this does not alter the nonzero structure.
1223 
1224    The user can set a value in the diagonal entry (or for the AIJ and
1225    row formats can optionally remove the main diagonal entry from the
1226    nonzero structure as well, by passing a null pointer as the final
1227    argument).
1228 
1229 .keywords: matrix, zero, rows, boundary conditions
1230 
1231 .seealso: MatZeroEntries(), MatGetSubMatrix(), MatGetSubMatrixInPlace()
1232 @*/
1233 int MatZeroRows(Mat mat,IS is, Scalar *diag)
1234 {
1235   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
1236   if (mat->ops.zerorows) return (*mat->ops.zerorows)(mat,is,diag);
1237   SETERRQ(PETSC_ERR_SUP,"MatZeroRows");
1238 }
1239 
1240 /*@
1241    MatGetSize - Returns the numbers of rows and columns in a matrix.
1242 
1243    Input Parameter:
1244 .  mat - the matrix
1245 
1246    Output Parameters:
1247 .  m - the number of global rows
1248 .  n - the number of global columns
1249 
1250 .keywords: matrix, dimension, size, rows, columns, global, get
1251 
1252 .seealso: MatGetLocalSize()
1253 @*/
1254 int MatGetSize(Mat mat,int *m,int* n)
1255 {
1256   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
1257   if (!m || !n) SETERRQ(1,"MatGetSize:Bad address for result");
1258   return (*mat->ops.getsize)(mat,m,n);
1259 }
1260 
1261 /*@
1262    MatGetLocalSize - Returns the number of rows and columns in a matrix
1263    stored locally.  This information may be implementation dependent, so
1264    use with care.
1265 
1266    Input Parameters:
1267 .  mat - the matrix
1268 
1269    Output Parameters:
1270 .  m - the number of local rows
1271 .  n - the number of local columns
1272 
1273 .keywords: matrix, dimension, size, local, rows, columns, get
1274 
1275 .seealso: MatGetSize()
1276 @*/
1277 int MatGetLocalSize(Mat mat,int *m,int* n)
1278 {
1279   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
1280   if (!m || !n) SETERRQ(1,"MatGetLocalSize:Bad address for result");
1281   return (*mat->ops.getlocalsize)(mat,m,n);
1282 }
1283 
1284 /*@
1285    MatGetOwnershipRange - Returns the range of matrix rows owned by
1286    this processor, assuming that the matrix is laid out with the first
1287    n1 rows on the first processor, the next n2 rows on the second, etc.
1288    For certain parallel layouts this range may not be well-defined.
1289 
1290    Input Parameters:
1291 .  mat - the matrix
1292 
1293    Output Parameters:
1294 .  m - the first local row
1295 .  n - one more then the last local row
1296 
1297 .keywords: matrix, get, range, ownership
1298 @*/
1299 int MatGetOwnershipRange(Mat mat,int *m,int* n)
1300 {
1301   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
1302   if (!m || !n) SETERRQ(1,"MatGetOwnershipRange:Bad address for result");
1303   if (mat->ops.getownershiprange) return (*mat->ops.getownershiprange)(mat,m,n);
1304   SETERRQ(PETSC_ERR_SUP,"MatGetOwnershipRange");
1305 }
1306 
1307 /*@
1308    MatILUFactorSymbolic - Performs symbolic ILU factorization of a matrix.
1309    Uses levels of fill only, not drop tolerance. Use MatLUFactorNumeric()
1310    to complete the factorization.
1311 
1312    Input Parameters:
1313 .  mat - the matrix
1314 .  row - row permutation
1315 .  column - column permutation
1316 .  fill - number of levels of fill
1317 .  f - expected fill as ratio of original fill
1318 
1319    Output Parameters:
1320 .  fact - puts factor
1321 
1322 .keywords: matrix, factor, incomplete, ILU, symbolic, fill
1323 
1324 .seealso: MatLUFactorSymbolic(), MatLUFactorNumeric()
1325 @*/
1326 int MatILUFactorSymbolic(Mat mat,IS row,IS col,double f,int fill,Mat *fact)
1327 {
1328   int ierr;
1329   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
1330   if (fill < 0) SETERRQ(1,"MatILUFactorSymbolic:Levels of fill negative");
1331   if (!fact) SETERRQ(1,"MatILUFactorSymbolic:Fact argument is missing");
1332   if (!mat->ops.ilufactorsymbolic) SETERRQ(PETSC_ERR_SUP,"MatILUFactorSymbolic");
1333   PLogEventBegin(MAT_ILUFactorSymbolic,mat,row,col,0);
1334   ierr = (*mat->ops.ilufactorsymbolic)(mat,row,col,f,fill,fact);
1335   CHKERRQ(ierr);
1336   PLogEventEnd(MAT_ILUFactorSymbolic,mat,row,col,0);
1337   return 0;
1338 }
1339 
1340 /*@
1341    MatIncompleteCholeskyFactorSymbolic - Performs symbolic incomplete
1342    Cholesky factorization for a symmetric matrix.  Use
1343    MatCholeskyFactorNumeric() to complete the factorization.
1344 
1345    Input Parameters:
1346 .  mat - the matrix
1347 .  perm - row and column permutation
1348 .  fill - levels of fill
1349 .  f - expected fill as ratio of original fill
1350 
1351    Output Parameter:
1352 .  fact - the factored matrix
1353 
1354    Note:  Currently only no-fill factorization is supported.
1355 
1356 .keywords: matrix, factor, incomplete, ICC, Cholesky, symbolic, fill
1357 
1358 .seealso: MatCholeskyFactorNumeric(), MatCholeskyFactor()
1359 @*/
1360 int MatIncompleteCholeskyFactorSymbolic(Mat mat,IS perm,double f,int fill,
1361                                         Mat *fact)
1362 {
1363   int ierr;
1364   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
1365   if (fill < 0) SETERRQ(1,"MatIncompleteCholeskyFactorSymbolic:Fill negative");
1366   if (!fact) SETERRQ(1,"MatIncompleteCholeskyFactorSymbolic:Missing fact argument");
1367   if (!mat->ops.incompletecholeskyfactorsymbolic)
1368     SETERRQ(PETSC_ERR_SUP,"MatIncompleteCholeskyFactorSymbolic");
1369   PLogEventBegin(MAT_IncompleteCholeskyFactorSymbolic,mat,perm,0,0);
1370   ierr = (*mat->ops.incompletecholeskyfactorsymbolic)(mat,perm,f,fill,fact);
1371   CHKERRQ(ierr);
1372   PLogEventEnd(MAT_IncompleteCholeskyFactorSymbolic,mat,perm,0,0);
1373   return 0;
1374 }
1375 
1376 /*@C
1377    MatGetArray - Returns a pointer to the element values in the matrix.
1378    This routine  is implementation dependent, and may not even work for
1379    certain matrix types.
1380 
1381    Input Parameter:
1382 .  mat - the matrix
1383 
1384    Output Parameter:
1385 .  v - the location of the values
1386 
1387 .keywords: matrix, array, elements, values
1388 @*/
1389 int MatGetArray(Mat mat,Scalar **v)
1390 {
1391   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
1392   if (!v) SETERRQ(1,"MatGetArray:Bad input, array pointer location");
1393   if (!mat->ops.getarray) SETERRQ(PETSC_ERR_SUP,"MatGetArraye");
1394   return (*mat->ops.getarray)(mat,v);
1395 }
1396 
1397 /*@C
1398    MatGetSubMatrix - Extracts a submatrix from a matrix. If submat points
1399                      to a valid matrix it may be reused.
1400 
1401    Input Parameters:
1402 .  mat - the matrix
1403 .  irow, icol - index sets of rows and columns to extract
1404 
1405    Output Parameter:
1406 .  submat - the submatrix
1407 
1408    Notes:
1409    MatGetSubMatrix() can be useful in setting boundary conditions.
1410 
1411 .keywords: matrix, get, submatrix, boundary conditions
1412 
1413 .seealso: MatZeroRows(), MatGetSubMatrixInPlace()
1414 @*/
1415 int MatGetSubMatrix(Mat mat,IS irow,IS icol,Mat *submat)
1416 {
1417   int ierr;
1418   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
1419   if (!mat->ops.getsubmatrix) SETERRQ(PETSC_ERR_SUP,"MatGetSubMatrix");
1420   PLogEventBegin(MAT_GetSubMatrix,mat,irow,icol,0);
1421   ierr = (*mat->ops.getsubmatrix)(mat,irow,icol,submat); CHKERRQ(ierr);
1422   PLogEventEnd(MAT_GetSubMatrix,mat,irow,icol,0);
1423   return 0;
1424 }
1425 
1426 /*@C
1427    MatGetSubMatrices - Extracts several submatrices from a matrix. If submat points
1428                      to an array of valid matrices it may be reused.
1429 
1430    Input Parameters:
1431 .  mat - the matrix
1432 .  irow, icol - index sets of rows and columns to extract
1433 
1434    Output Parameter:
1435 .  submat - the submatrices
1436 
1437 .keywords: matrix, get, submatrix
1438 
1439 .seealso: MatGetSubMatrix()
1440 @*/
1441 int MatGetSubMatrices(Mat mat,int n, IS *irow,IS *icol,Mat **submat)
1442 {
1443   int ierr;
1444   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
1445   if (!mat->ops.getsubmatrices) SETERRQ(PETSC_ERR_SUP,"MatGetSubMatrices");
1446   PLogEventBegin(MAT_GetSubMatrices,mat,0,0,0);
1447   ierr = (*mat->ops.getsubmatrices)(mat,n,irow,icol,submat); CHKERRQ(ierr);
1448   PLogEventEnd(MAT_GetSubMatrices,mat,0,0,0);
1449   return 0;
1450 }
1451 
1452 /*@
1453    MatGetSubMatrixInPlace - Extracts a submatrix from a matrix, returning
1454    the submatrix in place of the original matrix.
1455 
1456    Input Parameters:
1457 .  mat - the matrix
1458 .  irow, icol - index sets of rows and columns to extract
1459 
1460 .keywords: matrix, get, submatrix, boundary conditions, in-place
1461 
1462 .seealso: MatZeroRows(), MatGetSubMatrix()
1463 @*/
1464 int MatGetSubMatrixInPlace(Mat mat,IS irow,IS icol)
1465 {
1466   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
1467   if (!mat->ops.getsubmatrixinplace) SETERRQ(PETSC_ERR_SUP,"MatGetSubmatrixInPlace");
1468   return (*mat->ops.getsubmatrixinplace)(mat,irow,icol);
1469 }
1470 
1471 /*@
1472    MatGetType - Returns the type of the matrix, one of MATSEQDENSE, MATSEQAIJ, etc.
1473 
1474    Input Parameter:
1475 .  mat - the matrix
1476 
1477    Ouput Parameter:
1478 .  type - the matrix type
1479 
1480    Notes:
1481    The matrix types are given in petsc/include/mat.h
1482 
1483 .keywords: matrix, get, type
1484 
1485 .seealso:  MatGetName()
1486 @*/
1487 int MatGetType(Mat mat,MatType *type)
1488 {
1489   PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);
1490   *type = (MatType) mat->type;
1491   return 0;
1492 }
1493