xref: /petsc/src/mat/interface/matrix.c (revision bef22f13d3015174dc218807619adc0958744153)
1 #ifdef PETSC_RCS_HEADER
2 static char vcid[] = "$Id: matrix.c,v 1.287 1998/04/20 19:29:14 bsmith Exp curfman $";
3 #endif
4 
5 /*
6    This is where the abstract matrix operations are defined
7 */
8 
9 #include "src/mat/matimpl.h"        /*I "mat.h" I*/
10 #include "src/vec/vecimpl.h"
11 #include "pinclude/pviewer.h"
12 
13 #undef __FUNC__
14 #define __FUNC__ "MatGetRow"
15 /*@C
16    MatGetRow - Gets a row of a matrix.  You MUST call MatRestoreRow()
17    for each row that you get to ensure that your application does
18    not bleed memory.
19 
20    Not Collective
21 
22    Input Parameters:
23 +  mat - the matrix
24 -  row - the row to get
25 
26    Output Parameters:
27 +  ncols -  the number of nonzeros in the row
28 .  cols - if nonzero, the column numbers
29 -  vals - if nonzero, the values
30 
31    Notes:
32    This routine is provided for people who need to have direct access
33    to the structure of a matrix.  We hope that we provide enough
34    high-level matrix routines that few users will need it.
35 
36    MatGetRow() always returns 0-based column indices, regardless of
37    whether the internal representation is 0-based (default) or 1-based.
38 
39    For better efficiency, set cols and/or vals to PETSC_NULL if you do
40    not wish to extract these quantities.
41 
42    The user can only examine the values extracted with MatGetRow();
43    the values cannot be altered.  To change the matrix entries, one
44    must use MatSetValues().
45 
46    You can only have one call to MatGetRow() outstanding for a particular
47    matrix at a time.
48 
49    Fortran Notes:
50    The calling sequence from Fortran is
51 .vb
52    MatGetRow(matrix,row,ncols,cols,values,ierr)
53          Mat     matrix (input)
54          integer row    (input)
55          integer ncols  (output)
56          integer cols(maxcols) (output)
57          double precision (or double complex) values(maxcols) output
58 .ve
59    where maxcols >= maximum nonzeros in any row of the matrix.
60 
61    Caution:
62    Do not try to change the contents of the output arrays (cols and vals).
63    In some cases, this may corrupt the matrix.
64 
65 .keywords: matrix, row, get, extract
66 
67 .seealso: MatRestoreRow(), MatSetValues()
68 @*/
69 int MatGetRow(Mat mat,int row,int *ncols,int **cols,Scalar **vals)
70 {
71   int   ierr;
72 
73   PetscFunctionBegin;
74   PetscValidHeaderSpecific(mat,MAT_COOKIE);
75   PetscValidIntPointer(ncols);
76   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix");
77   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix");
78   if (!mat->ops->getrow) SETERRQ(PETSC_ERR_SUP,0,"");
79   PLogEventBegin(MAT_GetRow,mat,0,0,0);
80   ierr = (*mat->ops->getrow)(mat,row,ncols,cols,vals); CHKERRQ(ierr);
81   PLogEventEnd(MAT_GetRow,mat,0,0,0);
82   PetscFunctionReturn(0);
83 }
84 
85 #undef __FUNC__
86 #define __FUNC__ "MatRestoreRow"
87 /*@C
88    MatRestoreRow - Frees any temporary space allocated by MatGetRow().
89 
90    Not Collective
91 
92    Input Parameters:
93 +  mat - the matrix
94 .  row - the row to get
95 .  ncols, cols - the number of nonzeros and their columns
96 -  vals - if nonzero the column values
97 
98    Notes:
99    This routine should be called after you have finished examining the entries.
100 
101    Fortran Notes:
102    The calling sequence from Fortran is
103 .vb
104    MatRestoreRow(matrix,row,ncols,cols,values,ierr)
105       Mat     matrix (input)
106       integer row    (input)
107       integer ncols  (output)
108       integer cols(maxcols) (output)
109       double precision (or double complex) values(maxcols) output
110 .ve
111    Where maxcols >= maximum nonzeros in any row of the matrix.
112 
113    In Fortran MatRestoreRow() MUST be called after MatGetRow()
114    before another call to MatGetRow() can be made.
115 
116 .keywords: matrix, row, restore
117 
118 .seealso:  MatGetRow()
119 @*/
120 int MatRestoreRow(Mat mat,int row,int *ncols,int **cols,Scalar **vals)
121 {
122   int ierr;
123 
124   PetscFunctionBegin;
125   PetscValidHeaderSpecific(mat,MAT_COOKIE);
126   PetscValidIntPointer(ncols);
127   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix");
128   if (!mat->ops->restorerow) PetscFunctionReturn(0);
129   ierr = (*mat->ops->restorerow)(mat,row,ncols,cols,vals);CHKERRQ(ierr);
130   PetscFunctionReturn(0);
131 }
132 
133 #undef __FUNC__
134 #define __FUNC__ "MatView"
135 /*@C
136    MatView - Visualizes a matrix object.
137 
138    Input Parameters:
139 .  mat - the matrix
140 .  ptr - visualization context
141 
142    Collective on Mat unless Viewer is VIEWER_STDOUT_SELF
143 
144   Notes:
145    The available visualization contexts include
146 $     VIEWER_STDOUT_SELF - standard output (default)
147 $     VIEWER_STDOUT_WORLD - synchronized standard
148 $       output where only the first processor opens
149 $       the file.  All other processors send their
150 $       data to the first processor to print.
151 $     VIEWER_DRAWX_WORLD - graphical display of nonzero structure
152 
153    The user can open alternative vistualization contexts with
154 $    ViewerFileOpenASCII() - output to a specified file
155 $    ViewerFileOpenBinary() - output in binary to a
156 $         specified file; corresponding input uses MatLoad()
157 $    ViewerDrawOpenX() - output nonzero matrix structure to
158 $         an X window display
159 $    ViewerMatlabOpen() - output matrix to Matlab viewer.
160 $         Currently only the sequential dense and AIJ
161 $         matrix types support the Matlab viewer.
162 
163    The user can call ViewerSetFormat() to specify the output
164    format of ASCII printed objects (when using VIEWER_STDOUT_SELF,
165    VIEWER_STDOUT_WORLD and ViewerFileOpenASCII).  Available formats include
166 $    VIEWER_FORMAT_ASCII_DEFAULT - default, prints matrix contents
167 $    VIEWER_FORMAT_ASCII_MATLAB - Matlab format
168 $    VIEWER_FORMAT_ASCII_IMPL - implementation-specific format
169 $      (which is in many cases the same as the default)
170 $    VIEWER_FORMAT_ASCII_INFO - basic information about the matrix
171 $      size and structure (not the matrix entries)
172 $    VIEWER_FORMAT_ASCII_INFO_LONG - more detailed information about the
173 $      matrix structure
174 
175 .keywords: matrix, view, visualize, output, print, write, draw
176 
177 .seealso: ViewerSetFormat(), ViewerFileOpenASCII(), ViewerDrawOpenX(),
178           ViewerMatlabOpen(), ViewerFileOpenBinary(), MatLoad()
179 @*/
180 int MatView(Mat mat,Viewer viewer)
181 {
182   int          format, ierr, rows, cols;
183   FILE         *fd;
184   char         *cstr;
185   ViewerType   vtype;
186   MPI_Comm     comm = mat->comm;
187 
188   PetscFunctionBegin;
189   PetscValidHeaderSpecific(mat,MAT_COOKIE);
190   if (viewer) PetscValidHeaderSpecific(viewer,VIEWER_COOKIE);
191   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix");
192 
193   if (!viewer) {
194     viewer = VIEWER_STDOUT_SELF;
195   }
196 
197   ierr = ViewerGetType(viewer,&vtype);
198   if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
199     ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr);
200     ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
201     if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) {
202       PetscFPrintf(comm,fd,"Matrix Object:\n");
203       ierr = MatGetType(mat,PETSC_NULL,&cstr); CHKERRQ(ierr);
204       ierr = MatGetSize(mat,&rows,&cols); CHKERRQ(ierr);
205       PetscFPrintf(comm,fd,"  type=%s, rows=%d, cols=%d\n",cstr,rows,cols);
206       if (mat->ops->getinfo) {
207         MatInfo info;
208         ierr = MatGetInfo(mat,MAT_GLOBAL_SUM,&info); CHKERRQ(ierr);
209         PetscFPrintf(comm,fd,"  total: nonzeros=%d, allocated nonzeros=%d\n",
210                      (int)info.nz_used,(int)info.nz_allocated);
211       }
212     }
213   }
214   if (mat->ops->view) {ierr = (*mat->ops->view)(mat,viewer); CHKERRQ(ierr);}
215   PetscFunctionReturn(0);
216 }
217 
218 #undef __FUNC__
219 #define __FUNC__ "MatDestroy"
220 /*@C
221    MatDestroy - Frees space taken by a matrix.
222 
223    Input Parameter:
224 .  mat - the matrix
225 
226    Collective on Mat
227 
228 .keywords: matrix, destroy
229 @*/
230 int MatDestroy(Mat mat)
231 {
232   int ierr;
233 
234   PetscFunctionBegin;
235   PetscValidHeaderSpecific(mat,MAT_COOKIE);
236   if (--mat->refct > 0) PetscFunctionReturn(0);
237 
238   if (mat->mapping) {
239     ierr = ISLocalToGlobalMappingDestroy(mat->mapping); CHKERRQ(ierr);
240   }
241   if (mat->bmapping) {
242     ierr = ISLocalToGlobalMappingDestroy(mat->bmapping); CHKERRQ(ierr);
243   }
244   ierr = (*mat->ops->destroy)(mat); CHKERRQ(ierr);
245   PetscFunctionReturn(0);
246 }
247 
248 #undef __FUNC__
249 #define __FUNC__ "MatValid"
250 /*@
251    MatValid - Checks whether a matrix object is valid.
252 
253    Input Parameter:
254 .  m - the matrix to check
255 
256    Output Parameter:
257    flg - flag indicating matrix status, either
258 $     PETSC_TRUE if matrix is valid;
259 $     PETSC_FALSE otherwise.
260 
261    Collective on Mat
262 
263 .keywords: matrix, valid
264 @*/
265 int MatValid(Mat m,PetscTruth *flg)
266 {
267   PetscFunctionBegin;
268   PetscValidIntPointer(flg);
269   if (!m)                           *flg = PETSC_FALSE;
270   else if (m->cookie != MAT_COOKIE) *flg = PETSC_FALSE;
271   else                              *flg = PETSC_TRUE;
272   PetscFunctionReturn(0);
273 }
274 
275 #undef __FUNC__
276 #define __FUNC__ "MatSetValues"
277 /*@
278    MatSetValues - Inserts or adds a block of values into a matrix.
279    These values may be cached, so MatAssemblyBegin() and MatAssemblyEnd()
280    MUST be called after all calls to MatSetValues() have been completed.
281 
282    Input Parameters:
283 .  mat - the matrix
284 .  v - a logically two-dimensional array of values
285 .  m, idxm - the number of rows and their global indices
286 .  n, idxn - the number of columns and their global indices
287 .  addv - either ADD_VALUES or INSERT_VALUES, where
288 $     ADD_VALUES - adds values to any existing entries
289 $     INSERT_VALUES - replaces existing entries with new values
290 
291    Not Collective
292 
293    Notes:
294    By default the values, v, are row-oriented and unsorted.
295    See MatSetOption() for other options.
296 
297    Calls to MatSetValues() with the INSERT_VALUES and ADD_VALUES
298    options cannot be mixed without intervening calls to the assembly
299    routines.
300 
301    MatSetValues() uses 0-based row and column numbers in Fortran
302    as well as in C.
303 
304    Efficiency Alert:
305    The routine MatSetValuesBlocked() may offer much better efficiency
306    for users of block sparse formats (MATSEQBAIJ and MATMPIBAIJ).
307 
308 .keywords: matrix, insert, add, set, values
309 
310 .seealso: MatSetOption(), MatAssemblyBegin(), MatAssemblyEnd(), MatSetValuesBlocked()
311 @*/
312 int MatSetValues(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v,InsertMode addv)
313 {
314   int ierr;
315 
316   PetscFunctionBegin;
317   if (!m || !n) PetscFunctionReturn(0); /* no values to insert */
318   PetscValidHeaderSpecific(mat,MAT_COOKIE);
319   PetscValidIntPointer(idxm);
320   PetscValidIntPointer(idxn);
321   PetscValidScalarPointer(v);
322   if (mat->insertmode == NOT_SET_VALUES) {
323     mat->insertmode = addv;
324   }
325 #if defined(USE_PETSC_BOPT_g)
326   else if (mat->insertmode != addv) {
327     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,1,"Cannot mix add values and insert values");
328   }
329   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix");
330 #endif
331 
332   if (mat->assembled) {
333     mat->was_assembled = PETSC_TRUE;
334     mat->assembled     = PETSC_FALSE;
335   }
336   PLogEventBegin(MAT_SetValues,mat,0,0,0);
337   ierr = (*mat->ops->setvalues)(mat,m,idxm,n,idxn,v,addv);CHKERRQ(ierr);
338   PLogEventEnd(MAT_SetValues,mat,0,0,0);
339   PetscFunctionReturn(0);
340 }
341 
342 #undef __FUNC__
343 #define __FUNC__ "MatSetValuesBlocked"
344 /*@
345    MatSetValuesBlocked - Inserts or adds a block of values into a matrix.
346 
347    Input Parameters:
348 .  mat - the matrix
349 .  v - a logically two-dimensional array of values
350 .  m, idxm - the number of block rows and their global block indices
351 .  n, idxn - the number of block columns and their global block indices
352 .  addv - either ADD_VALUES or INSERT_VALUES, where
353 $     ADD_VALUES - adds values to any existing entries
354 $     INSERT_VALUES - replaces existing entries with new values
355 
356    Not Collective
357 
358    Notes:
359    By default the values, v, are row-oriented and unsorted. So the layout of
360    v is the same as for MatSetValues(). See MatSetOption() for other options.
361 
362    Calls to MatSetValuesBlocked() with the INSERT_VALUES and ADD_VALUES
363    options cannot be mixed without intervening calls to the assembly
364    routines.
365 
366    MatSetValuesBlocked() uses 0-based row and column numbers in Fortran
367    as well as in C.
368 
369    Each time an entry is set within a sparse matrix via MatSetValues(),
370    internal searching must be done to determine where to place the the
371    data in the matrix storage space.  By instead inserting blocks of
372    entries via MatSetValuesBlocked(), the overhead of matrix assembly is
373    reduced.
374 
375    Restrictions:
376    MatSetValuesBlocked() is currently supported only for the block AIJ
377    matrix format (MATSEQBAIJ and MATMPIBAIJ, which are created via
378    MatCreateSeqBAIJ() and MatCreateMPIBAIJ()).
379 
380 .keywords: matrix, insert, add, set, values
381 
382 .seealso: MatSetOption(), MatAssemblyBegin(), MatAssemblyEnd(), MatSetValues()
383 @*/
384 int MatSetValuesBlocked(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v,InsertMode addv)
385 {
386   int ierr;
387 
388   PetscFunctionBegin;
389   if (!m || !n) PetscFunctionReturn(0); /* no values to insert */
390   PetscValidHeaderSpecific(mat,MAT_COOKIE);
391   PetscValidIntPointer(idxm);
392   PetscValidIntPointer(idxn);
393   PetscValidScalarPointer(v);
394   if (mat->insertmode == NOT_SET_VALUES) {
395     mat->insertmode = addv;
396   }
397 #if defined(USE_PETSC_BOPT_g)
398   else if (mat->insertmode != addv) {
399     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,1,"Cannot mix add values and insert values");
400   }
401   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix");
402 #endif
403 
404   if (mat->assembled) {
405     mat->was_assembled = PETSC_TRUE;
406     mat->assembled     = PETSC_FALSE;
407   }
408   PLogEventBegin(MAT_SetValues,mat,0,0,0);
409   ierr = (*mat->ops->setvaluesblocked)(mat,m,idxm,n,idxn,v,addv);CHKERRQ(ierr);
410   PLogEventEnd(MAT_SetValues,mat,0,0,0);
411   PetscFunctionReturn(0);
412 }
413 
414 /*MC
415    MatSetValue - Set a single entry into a matrix.
416 
417    Input Parameters:
418 .  m - the matrix
419 .  row - the row location of the entry
420 .  col - the column location of the entry
421 .  value - the value to insert
422 .  mode - either INSERT_VALUES or ADD_VALUES
423 
424    Synopsis:
425    void MatSetValue(Mat m,int row,int col,Scalar value,InsertMode mode);
426 
427    Notes: For efficiency one should use MatSetValues() and set
428 several or many values simultaneously.
429 
430 .seealso: MatSetValues()
431 M*/
432 
433 #undef __FUNC__
434 #define __FUNC__ "MatGetValues"
435 /*@
436    MatGetValues - Gets a block of values from a matrix.
437 
438    Input Parameters:
439 .  mat - the matrix
440 .  v - a logically two-dimensional array for storing the values
441 .  m, idxm - the number of rows and their global indices
442 .  n, idxn - the number of columns and their global indices
443 
444    Not Collective; currently only returns a local block
445 
446    Notes:
447    The user must allocate space (m*n Scalars) for the values, v.
448    The values, v, are then returned in a row-oriented format,
449    analogous to that used by default in MatSetValues().
450 
451    MatGetValues() uses 0-based row and column numbers in
452    Fortran as well as in C.
453 
454    MatGetValues() requires that the matrix has been assembled
455    with MatAssemblyBegin()/MatAssemblyEnd().  Thus, calls to
456    MatSetValues() and MatGetValues() CANNOT be made in succession
457    without intermediate matrix assembly.
458 
459 .keywords: matrix, get, values
460 
461 .seealso: MatGetRow(), MatGetSubmatrices(), MatSetValues()
462 @*/
463 int MatGetValues(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v)
464 {
465   int ierr;
466 
467   PetscFunctionBegin;
468   PetscValidHeaderSpecific(mat,MAT_COOKIE);
469   PetscValidIntPointer(idxm);
470   PetscValidIntPointer(idxn);
471   PetscValidScalarPointer(v);
472   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix");
473   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix");
474   if (!mat->ops->getvalues) SETERRQ(PETSC_ERR_SUP,0,"");
475 
476   PLogEventBegin(MAT_GetValues,mat,0,0,0);
477   ierr = (*mat->ops->getvalues)(mat,m,idxm,n,idxn,v); CHKERRQ(ierr);
478   PLogEventEnd(MAT_GetValues,mat,0,0,0);
479   PetscFunctionReturn(0);
480 }
481 
482 #undef __FUNC__
483 #define __FUNC__ "MatSetLocalToGlobalMapping"
484 /*@
485    MatSetLocalToGlobalMapping - Sets a local-to-global numbering for use by
486    the routine MatSetValuesLocal() to allow users to insert matrix entries
487    using a local (per-processor) numbering.
488 
489    Input Parameters:
490 .  x - the matrix
491 .  mapping - mapping created with ISLocalToGlobalMappingCreate()
492              or ISLocalToGlobalMappingCreateIS()
493 
494    Not Collective
495 
496 .keywords: matrix, set, values, local, global, mapping
497 
498 .seealso:  MatAssemblyBegin(), MatAssemblyEnd(), MatSetValues(), MatSetValuesLocal()
499 @*/
500 int MatSetLocalToGlobalMapping(Mat x,ISLocalToGlobalMapping mapping)
501 {
502   PetscFunctionBegin;
503   PetscValidHeaderSpecific(x,MAT_COOKIE);
504   PetscValidHeaderSpecific(mapping,IS_LTOGM_COOKIE);
505 
506   if (x->mapping) {
507     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Mapping already set for matrix");
508   }
509 
510   x->mapping = mapping;
511   PetscObjectReference((PetscObject)mapping);
512   PetscFunctionReturn(0);
513 }
514 
515 #undef __FUNC__
516 #define __FUNC__ "MatSetLocalToGlobalMappingBlocked"
517 /*@
518    MatSetLocalToGlobalMappingBlocked - Sets a local-to-global numbering for use
519    by the routine MatSetValuesBlockedLocal() to allow users to insert matrix
520    entries using a local (per-processor) numbering.
521 
522    Input Parameters:
523 .  x - the matrix
524 .  mapping - mapping created with ISLocalToGlobalMappingCreate() or
525              ISLocalToGlobalMappingCreateIS()
526 
527    Not Collective
528 
529 .keywords: matrix, set, values, local ordering
530 
531 .seealso:  MatAssemblyBegin(), MatAssemblyEnd(), MatSetValues(), MatSetValuesBlockedLocal(),
532            MatSetValuesBlocked(), MatSetValuesLocal()
533 @*/
534 int MatSetLocalToGlobalMappingBlocked(Mat x,ISLocalToGlobalMapping mapping)
535 {
536   PetscFunctionBegin;
537   PetscValidHeaderSpecific(x,MAT_COOKIE);
538   PetscValidHeaderSpecific(mapping,IS_LTOGM_COOKIE);
539 
540   if (x->bmapping) {
541     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Mapping already set for matrix");
542   }
543 
544   x->bmapping = mapping;
545   PetscObjectReference((PetscObject)mapping);
546   PetscFunctionReturn(0);
547 }
548 
549 #undef __FUNC__
550 #define __FUNC__ "MatSetValuesLocal"
551 /*@
552    MatSetValuesLocal - Inserts or adds values into certain locations of a matrix,
553    using a local ordering of the nodes.
554 
555    Input Parameters:
556 .  x - the matrix
557 .  nrow, irow - number of rows and their local indices
558 .  ncol, icol - number of columns and their local indices
559 .  y -  a logically two-dimensional array of values
560 .  addv - either INSERT_VALUES or ADD_VALUES, where
561 $     ADD_VALUES - adds values to any existing entries
562 $     INSERT_VALUES - replaces existing entries with new values
563 
564    Not Collective
565 
566    Notes:
567    Before calling MatSetValuesLocal(), the user must first set the
568    local-to-global mapping by calling MatSetLocalToGlobalMapping().
569 
570    Calls to MatSetValuesLocal() with the INSERT_VALUES and ADD_VALUES
571    options cannot be mixed without intervening calls to the assembly
572    routines.
573 
574    These values may be cached, so MatAssemblyBegin() and MatAssemblyEnd()
575    MUST be called after all calls to MatSetValuesLocal() have been completed.
576 
577 .keywords: matrix, set, values, local ordering
578 
579 .seealso:  MatAssemblyBegin(), MatAssemblyEnd(), MatSetValues(), MatSetLocalToGlobalMapping()
580 @*/
581 int MatSetValuesLocal(Mat mat,int nrow,int *irow,int ncol, int *icol,Scalar *y,InsertMode addv)
582 {
583   int ierr,irowm[2048],icolm[2048];
584 
585   PetscFunctionBegin;
586   PetscValidHeaderSpecific(mat,MAT_COOKIE);
587   PetscValidIntPointer(irow);
588   PetscValidIntPointer(icol);
589   PetscValidScalarPointer(y);
590 
591   if (mat->insertmode == NOT_SET_VALUES) {
592     mat->insertmode = addv;
593   }
594 #if defined(USE_PETSC_BOPT_g)
595   else if (mat->insertmode != addv) {
596     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,1,"Cannot mix add values and insert values");
597   }
598   if (!mat->mapping) {
599     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Local to global never set with MatSetLocalToGlobalMapping");
600   }
601   if (nrow > 2048 || ncol > 2048) {
602     SETERRQ(PETSC_ERR_SUP,0,"Number column/row indices must be <= 2048");
603   }
604   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix");
605 #endif
606 
607   if (mat->assembled) {
608     mat->was_assembled = PETSC_TRUE;
609     mat->assembled     = PETSC_FALSE;
610   }
611   PLogEventBegin(MAT_SetValues,mat,0,0,0);
612   ierr = ISLocalToGlobalMappingApply(mat->mapping,nrow,irow,irowm); CHKERRQ(ierr);
613   ierr = ISLocalToGlobalMappingApply(mat->mapping,ncol,icol,icolm);CHKERRQ(ierr);
614   ierr = (*mat->ops->setvalues)(mat,nrow,irowm,ncol,icolm,y,addv);CHKERRQ(ierr);
615   PLogEventEnd(MAT_SetValues,mat,0,0,0);
616   PetscFunctionReturn(0);
617 }
618 
619 #undef __FUNC__
620 #define __FUNC__ "MatSetValuesBlockedLocal"
621 /*@
622    MatSetValuesBlockedLocal - Inserts or adds values into certain locations of a matrix,
623    using a local ordering of the nodes a block at a time.
624 
625    Input Parameters:
626 .  x - the matrix
627 .  nrow, irow - number of rows and their local indices
628 .  ncol, icol - number of columns and their local indices
629 .  y -  a logically two-dimensional array of values
630 .  addv - either INSERT_VALUES or ADD_VALUES, where
631 $     ADD_VALUES - adds values to any existing entries
632 $     INSERT_VALUES - replaces existing entries with new values
633 
634    Not Collective
635 
636    Notes:
637    Before calling MatSetValuesBlockedLocal(), the user must first set the
638    local-to-global mapping by calling MatSetLocalToGlobalMappingBlocked(),
639    where the mapping MUST be set for matrix blocks, not for matrix elements.
640 
641    Calls to MatSetValuesBlockedLocal() with the INSERT_VALUES and ADD_VALUES
642    options cannot be mixed without intervening calls to the assembly
643    routines.
644 
645    These values may be cached, so MatAssemblyBegin() and MatAssemblyEnd()
646    MUST be called after all calls to MatSetValuesBlockedLocal() have been completed.
647 
648 .keywords: matrix, set, values, blocked, local
649 
650 .seealso:  MatAssemblyBegin(), MatAssemblyEnd(), MatSetValuesLocal(), MatSetLocalToGlobalMappingBlocked()
651 @*/
652 int MatSetValuesBlockedLocal(Mat mat,int nrow,int *irow,int ncol,int *icol,Scalar *y,InsertMode addv)
653 {
654   int ierr,irowm[2048],icolm[2048];
655 
656   PetscFunctionBegin;
657   PetscValidHeaderSpecific(mat,MAT_COOKIE);
658   PetscValidIntPointer(irow);
659   PetscValidIntPointer(icol);
660   PetscValidScalarPointer(y);
661   if (mat->insertmode == NOT_SET_VALUES) {
662     mat->insertmode = addv;
663   }
664 #if defined(USE_PETSC_BOPT_g)
665   else if (mat->insertmode != addv) {
666     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,1,"Cannot mix add values and insert values");
667   }
668   if (!mat->bmapping) {
669     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Local to global never set with MatSetLocalToGlobalMappingBlocked");
670   }
671   if (nrow > 2048 || ncol > 2048) {
672     SETERRQ(PETSC_ERR_SUP,0,"Number column/row indices must be <= 2048");
673   }
674   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix");
675 #endif
676 
677   if (mat->assembled) {
678     mat->was_assembled = PETSC_TRUE;
679     mat->assembled     = PETSC_FALSE;
680   }
681   PLogEventBegin(MAT_SetValues,mat,0,0,0);
682   ierr = ISLocalToGlobalMappingApply(mat->bmapping,nrow,irow,irowm); CHKERRQ(ierr);
683   ierr = ISLocalToGlobalMappingApply(mat->bmapping,ncol,icol,icolm); CHKERRQ(ierr);
684   ierr = (*mat->ops->setvaluesblocked)(mat,nrow,irowm,ncol,icolm,y,addv);CHKERRQ(ierr);
685   PLogEventEnd(MAT_SetValues,mat,0,0,0);
686   PetscFunctionReturn(0);
687 }
688 
689 /* --------------------------------------------------------*/
690 #undef __FUNC__
691 #define __FUNC__ "MatMult"
692 /*@
693    MatMult - Computes the matrix-vector product, y = Ax.
694 
695    Input Parameters:
696 .  mat - the matrix
697 .  x   - the vector to be multilplied
698 
699    Output Parameters:
700 .  y - the result
701 
702    Collective on Mat and Vec
703 
704    Notes:
705    The vectors x and y cannot be the same.  I.e., one cannot
706    call MatMult(A,y,y).
707 
708 .keywords: matrix, multiply, matrix-vector product
709 
710 .seealso: MatMultTrans(), MatMultAdd(), MatMultTransAdd()
711 @*/
712 int MatMult(Mat mat,Vec x,Vec y)
713 {
714   int ierr;
715 
716   PetscFunctionBegin;
717   PetscValidHeaderSpecific(mat,MAT_COOKIE);
718   PetscValidHeaderSpecific(x,VEC_COOKIE);PetscValidHeaderSpecific(y,VEC_COOKIE);
719   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix");
720   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix");
721   if (x == y) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"x and y must be different vectors");
722   if (mat->N != x->N) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec x: global dim");
723   if (mat->M != y->N) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec y: global dim");
724   if (mat->m != y->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec y: local dim");
725 
726   PLogEventBegin(MAT_Mult,mat,x,y,0);
727   ierr = (*mat->ops->mult)(mat,x,y); CHKERRQ(ierr);
728   PLogEventEnd(MAT_Mult,mat,x,y,0);
729 
730   PetscFunctionReturn(0);
731 }
732 
733 #undef __FUNC__
734 #define __FUNC__ "MatMultTrans"
735 /*@
736    MatMultTrans - Computes matrix transpose times a vector.
737 
738    Input Parameters:
739 .  mat - the matrix
740 .  x   - the vector to be multilplied
741 
742    Output Parameters:
743 .  y - the result
744 
745    Collective on Mat and Vec
746 
747    Notes:
748    The vectors x and y cannot be the same.  I.e., one cannot
749    call MatMultTrans(A,y,y).
750 
751 .keywords: matrix, multiply, matrix-vector product, transpose
752 
753 .seealso: MatMult(), MatMultAdd(), MatMultTransAdd()
754 @*/
755 int MatMultTrans(Mat mat,Vec x,Vec y)
756 {
757   int ierr;
758 
759   PetscFunctionBegin;
760   PetscValidHeaderSpecific(mat,MAT_COOKIE);
761   PetscValidHeaderSpecific(x,VEC_COOKIE); PetscValidHeaderSpecific(y,VEC_COOKIE);
762   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix");
763   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix");
764   if (x == y) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"x and y must be different vectors");
765   if (mat->M != x->N) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec x: global dim");
766   if (mat->N != y->N) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec y: global dim");
767   PLogEventBegin(MAT_MultTrans,mat,x,y,0);
768   ierr = (*mat->ops->multtrans)(mat,x,y); CHKERRQ(ierr);
769   PLogEventEnd(MAT_MultTrans,mat,x,y,0);
770   PetscFunctionReturn(0);
771 }
772 
773 #undef __FUNC__
774 #define __FUNC__ "MatMultAdd"
775 /*@
776     MatMultAdd -  Computes v3 = v2 + A * v1.
777 
778     Input Parameters:
779 .   mat - the matrix
780 .   v1, v2 - the vectors
781 
782     Output Parameters:
783 .   v3 - the result
784 
785    Collective on Mat and Vec
786 
787    Notes:
788    The vectors v1 and v3 cannot be the same.  I.e., one cannot
789    call MatMultAdd(A,v1,v2,v1).
790 
791 .keywords: matrix, multiply, matrix-vector product, add
792 
793 .seealso: MatMultTrans(), MatMult(), MatMultTransAdd()
794 @*/
795 int MatMultAdd(Mat mat,Vec v1,Vec v2,Vec v3)
796 {
797   int ierr;
798 
799   PetscFunctionBegin;
800   PetscValidHeaderSpecific(mat,MAT_COOKIE);PetscValidHeaderSpecific(v1,VEC_COOKIE);
801   PetscValidHeaderSpecific(v2,VEC_COOKIE); PetscValidHeaderSpecific(v3,VEC_COOKIE);
802   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix");
803   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix");
804   if (mat->N != v1->N) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec v1: global dim");
805   if (mat->M != v2->N) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec v2: global dim");
806   if (mat->M != v3->N) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec v3: global dim");
807   if (mat->m != v3->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec v3: local dim");
808   if (mat->m != v2->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec v2: local dim");
809   if (v1 == v3) SETERRQ(PETSC_ERR_ARG_IDN,0,"v1 and v3 must be different vectors");
810 
811   PLogEventBegin(MAT_MultAdd,mat,v1,v2,v3);
812   ierr = (*mat->ops->multadd)(mat,v1,v2,v3); CHKERRQ(ierr);
813   PLogEventEnd(MAT_MultAdd,mat,v1,v2,v3);
814   PetscFunctionReturn(0);
815 }
816 
817 #undef __FUNC__
818 #define __FUNC__ "MatMultTransAdd"
819 /*@
820    MatMultTransAdd - Computes v3 = v2 + A' * v1.
821 
822    Input Parameters:
823 .  mat - the matrix
824 .  v1, v2 - the vectors
825 
826    Output Parameters:
827 .  v3 - the result
828 
829    Collective on Mat and Vec
830 
831    Notes:
832    The vectors v1 and v3 cannot be the same.  I.e., one cannot
833    call MatMultTransAdd(A,v1,v2,v1).
834 
835 .keywords: matrix, multiply, matrix-vector product, transpose, add
836 
837 .seealso: MatMultTrans(), MatMultAdd(), MatMult()
838 @*/
839 int MatMultTransAdd(Mat mat,Vec v1,Vec v2,Vec v3)
840 {
841   int ierr;
842 
843   PetscFunctionBegin;
844   PetscValidHeaderSpecific(mat,MAT_COOKIE);PetscValidHeaderSpecific(v1,VEC_COOKIE);
845   PetscValidHeaderSpecific(v2,VEC_COOKIE);PetscValidHeaderSpecific(v3,VEC_COOKIE);
846   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix");
847   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix");
848   if (!mat->ops->multtransadd) SETERRQ(PETSC_ERR_SUP,0,"");
849   if (v1 == v3) SETERRQ(PETSC_ERR_ARG_IDN,0,"v1 and v3 must be different vectors");
850   if (mat->M != v1->N) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec v1: global dim");
851   if (mat->N != v2->N) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec v2: global dim");
852   if (mat->N != v3->N) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec v3: global dim");
853 
854   PLogEventBegin(MAT_MultTransAdd,mat,v1,v2,v3);
855   ierr = (*mat->ops->multtransadd)(mat,v1,v2,v3); CHKERRQ(ierr);
856   PLogEventEnd(MAT_MultTransAdd,mat,v1,v2,v3);
857   PetscFunctionReturn(0);
858 }
859 /* ------------------------------------------------------------*/
860 #undef __FUNC__
861 #define __FUNC__ "MatGetInfo"
862 /*@C
863    MatGetInfo - Returns information about matrix storage (number of
864    nonzeros, memory, etc.).
865 
866    Input Parameters:
867 .  mat - the matrix
868 
869    Output Parameters:
870 .  flag - flag indicating the type of parameters to be returned
871 $    flag = MAT_LOCAL: local matrix
872 $    flag = MAT_GLOBAL_MAX: maximum over all processors
873 $    flag = MAT_GLOBAL_SUM: sum over all processors
874 .  info - matrix information context
875 
876     Collective on Mat if MAT_GLOBAL_MAX or MAT_GLOBAL_SUM is passed in
877 
878    Notes:
879    The MatInfo context contains a variety of matrix data, including
880    number of nonzeros allocated and used, number of mallocs during
881    matrix assembly, etc.  Additional information for factored matrices
882    is provided (such as the fill ratio, number of mallocs during
883    factorization, etc.).  Much of this info is printed to STDOUT
884    when using the runtime options
885 $       -log_info -mat_view_info
886 
887    Example for C/C++ Users:
888    See the file ${PETSC_DIR}/include/mat.h for a complete list of
889    data within the MatInfo context.  For example,
890 $
891 $      MatInfo info;
892 $      Mat     A;
893 $      double  mal, nz_a, nz_u;
894 $
895 $      MatGetInfo(A,MAT_LOCAL,&info);
896 $      mal  = info.mallocs;
897 $      nz_a = info.nz_allocated;
898 $
899 
900    Example for Fortran Users:
901    Fortran users should declare info as a double precision
902    array of dimension MAT_INFO_SIZE, and then extract the parameters
903    of interest.  See the file ${PETSC_DIR}/include/finclude/mat.h
904    a complete list of parameter names.
905 $
906 $      double  precision info(MAT_INFO_SIZE)
907 $      double  precision mal, nz_a
908 $      Mat     A
909 $      integer ierr
910 $
911 $      call MatGetInfo(A,MAT_LOCAL,info,ierr)
912 $      mal = info(MAT_INFO_MALLOCS)
913 $      nz_a = info(MAT_INFO_NZ_ALLOCATED)
914 $
915 
916 .keywords: matrix, get, info, storage, nonzeros, memory, fill
917 @*/
918 int MatGetInfo(Mat mat,MatInfoType flag,MatInfo *info)
919 {
920   int ierr;
921 
922   PetscFunctionBegin;
923   PetscValidHeaderSpecific(mat,MAT_COOKIE);
924   PetscValidPointer(info);
925   if (!mat->ops->getinfo) SETERRQ(PETSC_ERR_SUP,0,"");
926   ierr = (*mat->ops->getinfo)(mat,flag,info);CHKERRQ(ierr);
927   PetscFunctionReturn(0);
928 }
929 
930 /* ----------------------------------------------------------*/
931 #undef __FUNC__
932 #define __FUNC__ "MatILUDTFactor"
933 /*@
934    MatILUDTFactor - Performs a drop tolerance ILU factorization.
935 
936    Input Parameters:
937 .  mat - the matrix
938 .  dt  - the drop tolerance
939 .  maxnz - the maximum number of nonzeros per row allowed?
940 .  row - row permutation
941 .  col - column permutation
942 
943    Output Parameters:
944 .  fact - the factored matrix
945 
946    Collective on Mat
947 
948 .keywords: matrix, factor, LU, in-place
949 
950 .seealso: MatLUFactorSymbolic(), MatLUFactorNumeric(), MatCholeskyFactor()
951 @*/
952 int MatILUDTFactor(Mat mat,double dt,int maxnz,IS row,IS col,Mat *fact)
953 {
954   int ierr;
955 
956   PetscFunctionBegin;
957   PetscValidHeaderSpecific(mat,MAT_COOKIE);
958   PetscValidPointer(fact);
959   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix");
960   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix");
961   if (!mat->ops->iludtfactor) SETERRQ(PETSC_ERR_SUP,0,"");
962 
963   PLogEventBegin(MAT_ILUFactor,mat,row,col,0);
964   ierr = (*mat->ops->iludtfactor)(mat,dt,maxnz,row,col,fact); CHKERRQ(ierr);
965   PLogEventEnd(MAT_ILUFactor,mat,row,col,0);
966 
967   PetscFunctionReturn(0);
968 }
969 
970 #undef __FUNC__
971 #define __FUNC__ "MatLUFactor"
972 /*@
973    MatLUFactor - Performs in-place LU factorization of matrix.
974 
975    Input Parameters:
976 .  mat - the matrix
977 .  row - row permutation
978 .  col - column permutation
979 .  f - expected fill as ratio of original fill.
980 
981    Collective on Mat
982 
983 .keywords: matrix, factor, LU, in-place
984 
985 .seealso: MatLUFactorSymbolic(), MatLUFactorNumeric(), MatCholeskyFactor()
986 @*/
987 int MatLUFactor(Mat mat,IS row,IS col,double f)
988 {
989   int ierr;
990 
991   PetscFunctionBegin;
992   PetscValidHeaderSpecific(mat,MAT_COOKIE);
993   if (mat->M != mat->N) SETERRQ(PETSC_ERR_ARG_WRONG,0,"matrix must be square");
994   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix");
995   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix");
996   if (!mat->ops->lufactor) SETERRQ(PETSC_ERR_SUP,0,"");
997 
998   PLogEventBegin(MAT_LUFactor,mat,row,col,0);
999   ierr = (*mat->ops->lufactor)(mat,row,col,f); CHKERRQ(ierr);
1000   PLogEventEnd(MAT_LUFactor,mat,row,col,0);
1001   PetscFunctionReturn(0);
1002 }
1003 
1004 #undef __FUNC__
1005 #define __FUNC__ "MatILUFactor"
1006 /*@
1007    MatILUFactor - Performs in-place ILU factorization of matrix.
1008 
1009    Input Parameters:
1010 .  mat - the matrix
1011 .  row - row permutation
1012 .  col - column permutation
1013 .  f - expected fill as ratio of original fill.
1014 .  level - number of levels of fill.
1015 
1016    Collective on Mat
1017 
1018    Notes:
1019    Probably really in-place only when level of fill is zero, otherwise allocates
1020    new space to store factored matrix and deletes previous memory.
1021 
1022 .keywords: matrix, factor, ILU, in-place
1023 
1024 .seealso: MatILUFactorSymbolic(), MatLUFactorNumeric(), MatCholeskyFactor()
1025 @*/
1026 int MatILUFactor(Mat mat,IS row,IS col,double f,int level)
1027 {
1028   int ierr;
1029 
1030   PetscFunctionBegin;
1031   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1032   if (mat->M != mat->N) SETERRQ(PETSC_ERR_ARG_WRONG,0,"matrix must be square");
1033   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix");
1034   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix");
1035   if (!mat->ops->ilufactor) SETERRQ(PETSC_ERR_SUP,0,"");
1036 
1037   PLogEventBegin(MAT_ILUFactor,mat,row,col,0);
1038   ierr = (*mat->ops->ilufactor)(mat,row,col,f,level); CHKERRQ(ierr);
1039   PLogEventEnd(MAT_ILUFactor,mat,row,col,0);
1040   PetscFunctionReturn(0);
1041 }
1042 
1043 #undef __FUNC__
1044 #define __FUNC__ "MatLUFactorSymbolic"
1045 /*@
1046    MatLUFactorSymbolic - Performs symbolic LU factorization of matrix.
1047    Call this routine before calling MatLUFactorNumeric().
1048 
1049    Input Parameters:
1050 .  mat - the matrix
1051 .  row, col - row and column permutations
1052 .  f - expected fill as ratio of the original number of nonzeros,
1053        for example 3.0; choosing this parameter well can result in
1054        more efficient use of time and space. Run with the option -log_info
1055        to determine an optimal value to use
1056 
1057    Output Parameter:
1058 .  fact - new matrix that has been symbolically factored
1059 
1060    Collective on Mat
1061 
1062    Notes:
1063    See the file ${PETSC_DIR}/Performance for additional information about
1064    choosing the fill factor for better efficiency.
1065 
1066 .keywords: matrix, factor, LU, symbolic, fill
1067 
1068 .seealso: MatLUFactor(), MatLUFactorNumeric(), MatCholeskyFactor()
1069 @*/
1070 int MatLUFactorSymbolic(Mat mat,IS row,IS col,double f,Mat *fact)
1071 {
1072   int ierr;
1073 
1074   PetscFunctionBegin;
1075   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1076   if (mat->M != mat->N) SETERRQ(PETSC_ERR_ARG_WRONG,0,"matrix must be square");
1077   PetscValidPointer(fact);
1078   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix");
1079   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix");
1080   if (!mat->ops->lufactorsymbolic) SETERRQ(PETSC_ERR_SUP,0,"");
1081 
1082   PLogEventBegin(MAT_LUFactorSymbolic,mat,row,col,0);
1083   ierr = (*mat->ops->lufactorsymbolic)(mat,row,col,f,fact); CHKERRQ(ierr);
1084   PLogEventEnd(MAT_LUFactorSymbolic,mat,row,col,0);
1085   PetscFunctionReturn(0);
1086 }
1087 
1088 #undef __FUNC__
1089 #define __FUNC__ "MatLUFactorNumeric"
1090 /*@
1091    MatLUFactorNumeric - Performs numeric LU factorization of a matrix.
1092    Call this routine after first calling MatLUFactorSymbolic().
1093 
1094    Input Parameters:
1095 .  mat - the matrix
1096 .  row, col - row and  column permutations
1097 
1098    Output Parameters:
1099 .  fact - symbolically factored matrix that must have been generated
1100           by MatLUFactorSymbolic()
1101 
1102    Collective on Mat
1103 
1104    Notes:
1105    See MatLUFactor() for in-place factorization.  See
1106    MatCholeskyFactorNumeric() for the symmetric, positive definite case.
1107 
1108 .keywords: matrix, factor, LU, numeric
1109 
1110 .seealso: MatLUFactorSymbolic(), MatLUFactor(), MatCholeskyFactor()
1111 @*/
1112 int MatLUFactorNumeric(Mat mat,Mat *fact)
1113 {
1114   int ierr,flg;
1115 
1116   PetscFunctionBegin;
1117   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1118   PetscValidPointer(fact);
1119   PetscValidHeaderSpecific(*fact,MAT_COOKIE);
1120   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix");
1121   if (mat->M != (*fact)->M || mat->N != (*fact)->N) {
1122     SETERRQ(PETSC_ERR_ARG_SIZ,0,"Mat mat,Mat *fact: global dimensions are different");
1123   }
1124   if (!(*fact)->ops->lufactornumeric) SETERRQ(PETSC_ERR_SUP,0,"");
1125 
1126   PLogEventBegin(MAT_LUFactorNumeric,mat,*fact,0,0);
1127   ierr = (*(*fact)->ops->lufactornumeric)(mat,fact); CHKERRQ(ierr);
1128   PLogEventEnd(MAT_LUFactorNumeric,mat,*fact,0,0);
1129   ierr = OptionsHasName(PETSC_NULL,"-mat_view_draw",&flg); CHKERRQ(ierr);
1130   if (flg) {
1131     ierr = OptionsHasName(0,"-mat_view_contour",&flg); CHKERRQ(ierr);
1132     if (flg) {
1133       ViewerPushFormat(VIEWER_DRAWX_(mat->comm),VIEWER_FORMAT_DRAW_CONTOUR,0);CHKERRQ(ierr);
1134     }
1135     ierr = MatView(*fact,VIEWER_DRAWX_(mat->comm)); CHKERRQ(ierr);
1136     ierr = ViewerFlush(VIEWER_DRAWX_(mat->comm)); CHKERRQ(ierr);
1137     if (flg) {
1138       ViewerPopFormat(VIEWER_DRAWX_(mat->comm));CHKERRQ(ierr);
1139     }
1140   }
1141   PetscFunctionReturn(0);
1142 }
1143 
1144 #undef __FUNC__
1145 #define __FUNC__ "MatCholeskyFactor"
1146 /*@
1147    MatCholeskyFactor - Performs in-place Cholesky factorization of a
1148    symmetric matrix.
1149 
1150    Input Parameters:
1151 .  mat - the matrix
1152 .  perm - row and column permutations
1153 .  f - expected fill as ratio of original fill
1154 
1155    Collective on Mat
1156 
1157    Notes:
1158    See MatLUFactor() for the nonsymmetric case.  See also
1159    MatCholeskyFactorSymbolic(), and MatCholeskyFactorNumeric().
1160 
1161 .keywords: matrix, factor, in-place, Cholesky
1162 
1163 .seealso: MatLUFactor(), MatCholeskyFactorSymbolic(), MatCholeskyFactorNumeric()
1164 @*/
1165 int MatCholeskyFactor(Mat mat,IS perm,double f)
1166 {
1167   int ierr;
1168 
1169   PetscFunctionBegin;
1170   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1171   if (mat->M != mat->N) SETERRQ(PETSC_ERR_ARG_WRONG,0,"matrix must be square");
1172   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix");
1173   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix");
1174   if (!mat->ops->choleskyfactor) SETERRQ(PETSC_ERR_SUP,0,"");
1175 
1176   PLogEventBegin(MAT_CholeskyFactor,mat,perm,0,0);
1177   ierr = (*mat->ops->choleskyfactor)(mat,perm,f); CHKERRQ(ierr);
1178   PLogEventEnd(MAT_CholeskyFactor,mat,perm,0,0);
1179   PetscFunctionReturn(0);
1180 }
1181 
1182 #undef __FUNC__
1183 #define __FUNC__ "MatCholeskyFactorSymbolic"
1184 /*@
1185    MatCholeskyFactorSymbolic - Performs symbolic Cholesky factorization
1186    of a symmetric matrix.
1187 
1188    Input Parameters:
1189 .  mat - the matrix
1190 .  perm - row and column permutations
1191 .  f - expected fill as ratio of original
1192 
1193    Output Parameter:
1194 .  fact - the factored matrix
1195 
1196    Collective on Mat
1197 
1198    Notes:
1199    See MatLUFactorSymbolic() for the nonsymmetric case.  See also
1200    MatCholeskyFactor() and MatCholeskyFactorNumeric().
1201 
1202 .keywords: matrix, factor, factorization, symbolic, Cholesky
1203 
1204 .seealso: MatLUFactorSymbolic(), MatCholeskyFactor(), MatCholeskyFactorNumeric()
1205 @*/
1206 int MatCholeskyFactorSymbolic(Mat mat,IS perm,double f,Mat *fact)
1207 {
1208   int ierr;
1209 
1210   PetscFunctionBegin;
1211   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1212   PetscValidPointer(fact);
1213   if (mat->M != mat->N) SETERRQ(PETSC_ERR_ARG_WRONG,0,"matrix must be square");
1214   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix");
1215   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix");
1216   if (!mat->ops->choleskyfactorsymbolic) SETERRQ(PETSC_ERR_SUP,0,"");
1217 
1218   PLogEventBegin(MAT_CholeskyFactorSymbolic,mat,perm,0,0);
1219   ierr = (*mat->ops->choleskyfactorsymbolic)(mat,perm,f,fact); CHKERRQ(ierr);
1220   PLogEventEnd(MAT_CholeskyFactorSymbolic,mat,perm,0,0);
1221   PetscFunctionReturn(0);
1222 }
1223 
1224 #undef __FUNC__
1225 #define __FUNC__ "MatCholeskyFactorNumeric"
1226 /*@
1227    MatCholeskyFactorNumeric - Performs numeric Cholesky factorization
1228    of a symmetric matrix. Call this routine after first calling
1229    MatCholeskyFactorSymbolic().
1230 
1231    Input Parameter:
1232 .  mat - the initial matrix
1233 
1234    Output Parameter:
1235 .  fact - the factored matrix
1236 
1237    Collective on Mat
1238 
1239 .keywords: matrix, factor, numeric, Cholesky
1240 
1241 .seealso: MatCholeskyFactorSymbolic(), MatCholeskyFactor(), MatLUFactorNumeric()
1242 @*/
1243 int MatCholeskyFactorNumeric(Mat mat,Mat *fact)
1244 {
1245   int ierr;
1246 
1247   PetscFunctionBegin;
1248   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1249   PetscValidPointer(fact);
1250   if (!mat->ops->choleskyfactornumeric) SETERRQ(PETSC_ERR_SUP,0,"");
1251   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix");
1252   if (mat->M != (*fact)->M || mat->N != (*fact)->N)
1253     SETERRQ(PETSC_ERR_ARG_SIZ,0,"Mat mat,Mat *fact: global dim");
1254 
1255   PLogEventBegin(MAT_CholeskyFactorNumeric,mat,*fact,0,0);
1256   ierr = (*mat->ops->choleskyfactornumeric)(mat,fact); CHKERRQ(ierr);
1257   PLogEventEnd(MAT_CholeskyFactorNumeric,mat,*fact,0,0);
1258   PetscFunctionReturn(0);
1259 }
1260 
1261 /* ----------------------------------------------------------------*/
1262 #undef __FUNC__
1263 #define __FUNC__ "MatSolve"
1264 /*@
1265    MatSolve - Solves A x = b, given a factored matrix.
1266 
1267    Input Parameters:
1268 .  mat - the factored matrix
1269 .  b - the right-hand-side vector
1270 
1271    Output Parameter:
1272 .  x - the result vector
1273 
1274    Collective on Mat and Vec
1275 
1276    Notes:
1277    The vectors b and x cannot be the same.  I.e., one cannot
1278    call MatSolve(A,x,x).
1279 
1280 .keywords: matrix, linear system, solve, LU, Cholesky, triangular solve
1281 
1282 .seealso: MatSolveAdd(), MatSolveTrans(), MatSolveTransAdd()
1283 @*/
1284 int MatSolve(Mat mat,Vec b,Vec x)
1285 {
1286   int ierr;
1287 
1288   PetscFunctionBegin;
1289   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1290   PetscValidHeaderSpecific(b,VEC_COOKIE); PetscValidHeaderSpecific(x,VEC_COOKIE);
1291   if (x == b) SETERRQ(PETSC_ERR_ARG_IDN,0,"x and b must be different vectors");
1292   if (!mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Unfactored matrix");
1293   if (mat->N != x->N) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec x: global dim");
1294   if (mat->M != b->N) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec b: global dim");
1295   if (mat->m != b->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec b: local dim");
1296   if (mat->M == 0 && mat->N == 0) PetscFunctionReturn(0);
1297 
1298   if (!mat->ops->solve) SETERRQ(PETSC_ERR_SUP,0,"");
1299   PLogEventBegin(MAT_Solve,mat,b,x,0);
1300   ierr = (*mat->ops->solve)(mat,b,x); CHKERRQ(ierr);
1301   PLogEventEnd(MAT_Solve,mat,b,x,0);
1302   PetscFunctionReturn(0);
1303 }
1304 
1305 #undef __FUNC__
1306 #define __FUNC__ "MatForwardSolve"
1307 /* @
1308    MatForwardSolve - Solves L x = b, given a factored matrix, A = LU.
1309 
1310    Input Parameters:
1311 .  mat - the factored matrix
1312 .  b - the right-hand-side vector
1313 
1314    Output Parameter:
1315 .  x - the result vector
1316 
1317    Notes:
1318    MatSolve() should be used for most applications, as it performs
1319    a forward solve followed by a backward solve.
1320 
1321    The vectors b and x cannot be the same.  I.e., one cannot
1322    call MatForwardSolve(A,x,x).
1323 
1324 .keywords: matrix, forward, LU, Cholesky, triangular solve
1325 
1326 .seealso: MatSolve(), MatBackwardSolve()
1327 @ */
1328 int MatForwardSolve(Mat mat,Vec b,Vec x)
1329 {
1330   int ierr;
1331 
1332   PetscFunctionBegin;
1333   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1334   PetscValidHeaderSpecific(b,VEC_COOKIE);  PetscValidHeaderSpecific(x,VEC_COOKIE);
1335   if (x == b) SETERRQ(PETSC_ERR_ARG_IDN,0,"x and b must be different vectors");
1336   if (!mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Unfactored matrix");
1337   if (!mat->ops->forwardsolve) SETERRQ(PETSC_ERR_SUP,0,"");
1338   if (mat->N != x->N) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec x: global dim");
1339   if (mat->M != b->N) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec b: global dim");
1340   if (mat->m != b->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec b: local dim");
1341 
1342   PLogEventBegin(MAT_ForwardSolve,mat,b,x,0);
1343   ierr = (*mat->ops->forwardsolve)(mat,b,x); CHKERRQ(ierr);
1344   PLogEventEnd(MAT_ForwardSolve,mat,b,x,0);
1345   PetscFunctionReturn(0);
1346 }
1347 
1348 #undef __FUNC__
1349 #define __FUNC__ "MatBackwardSolve"
1350 /* @
1351    MatBackwardSolve - Solves U x = b, given a factored matrix, A = LU.
1352 
1353    Input Parameters:
1354 .  mat - the factored matrix
1355 .  b - the right-hand-side vector
1356 
1357    Output Parameter:
1358 .  x - the result vector
1359 
1360    Notes:
1361    MatSolve() should be used for most applications, as it performs
1362    a forward solve followed by a backward solve.
1363 
1364    The vectors b and x cannot be the same.  I.e., one cannot
1365    call MatBackwardSolve(A,x,x).
1366 
1367 .keywords: matrix, backward, LU, Cholesky, triangular solve
1368 
1369 .seealso: MatSolve(), MatForwardSolve()
1370 @ */
1371 int MatBackwardSolve(Mat mat,Vec b,Vec x)
1372 {
1373   int ierr;
1374 
1375   PetscFunctionBegin;
1376   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1377   PetscValidHeaderSpecific(b,VEC_COOKIE);  PetscValidHeaderSpecific(x,VEC_COOKIE);
1378   if (x == b) SETERRQ(PETSC_ERR_ARG_IDN,0,"x and b must be different vectors");
1379   if (!mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Unfactored matrix");
1380   if (!mat->ops->backwardsolve) SETERRQ(PETSC_ERR_SUP,0,"");
1381   if (mat->N != x->N) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec x: global dim");
1382   if (mat->M != b->N) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec b: global dim");
1383   if (mat->m != b->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec b: local dim");
1384 
1385   PLogEventBegin(MAT_BackwardSolve,mat,b,x,0);
1386   ierr = (*mat->ops->backwardsolve)(mat,b,x); CHKERRQ(ierr);
1387   PLogEventEnd(MAT_BackwardSolve,mat,b,x,0);
1388   PetscFunctionReturn(0);
1389 }
1390 
1391 #undef __FUNC__
1392 #define __FUNC__ "MatSolveAdd"
1393 /*@
1394    MatSolveAdd - Computes x = y + inv(A)*b, given a factored matrix.
1395 
1396    Input Parameters:
1397 .  mat - the factored matrix
1398 .  b - the right-hand-side vector
1399 .  y - the vector to be added to
1400 
1401    Output Parameter:
1402 .  x - the result vector
1403 
1404    Collective on Mat and Vec
1405 
1406    Notes:
1407    The vectors b and x cannot be the same.  I.e., one cannot
1408    call MatSolveAdd(A,x,y,x).
1409 
1410 .keywords: matrix, linear system, solve, LU, Cholesky, add
1411 
1412 .seealso: MatSolve(), MatSolveTrans(), MatSolveTransAdd()
1413 @*/
1414 int MatSolveAdd(Mat mat,Vec b,Vec y,Vec x)
1415 {
1416   Scalar one = 1.0;
1417   Vec    tmp;
1418   int    ierr;
1419 
1420   PetscFunctionBegin;
1421   PetscValidHeaderSpecific(mat,MAT_COOKIE);PetscValidHeaderSpecific(y,VEC_COOKIE);
1422   PetscValidHeaderSpecific(b,VEC_COOKIE);  PetscValidHeaderSpecific(x,VEC_COOKIE);
1423   if (x == b) SETERRQ(PETSC_ERR_ARG_IDN,0,"x and b must be different vectors");
1424   if (!mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Unfactored matrix");
1425   if (mat->N != x->N) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec x: global dim");
1426   if (mat->M != b->N) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec b: global dim");
1427   if (mat->M != y->N) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec y: global dim");
1428   if (mat->m != b->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec b: local dim");
1429   if (x->n != y->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Vec x,Vec y: local dim");
1430 
1431   PLogEventBegin(MAT_SolveAdd,mat,b,x,y);
1432   if (mat->ops->solveadd)  {
1433     ierr = (*mat->ops->solveadd)(mat,b,y,x); CHKERRQ(ierr);
1434   }
1435   else {
1436     /* do the solve then the add manually */
1437     if (x != y) {
1438       ierr = MatSolve(mat,b,x); CHKERRQ(ierr);
1439       ierr = VecAXPY(&one,y,x); CHKERRQ(ierr);
1440     } else {
1441       ierr = VecDuplicate(x,&tmp); CHKERRQ(ierr);
1442       PLogObjectParent(mat,tmp);
1443       ierr = VecCopy(x,tmp); CHKERRQ(ierr);
1444       ierr = MatSolve(mat,b,x); CHKERRQ(ierr);
1445       ierr = VecAXPY(&one,tmp,x); CHKERRQ(ierr);
1446       ierr = VecDestroy(tmp); CHKERRQ(ierr);
1447     }
1448   }
1449   PLogEventEnd(MAT_SolveAdd,mat,b,x,y);
1450   PetscFunctionReturn(0);
1451 }
1452 
1453 #undef __FUNC__
1454 #define __FUNC__ "MatSolveTrans"
1455 /*@
1456    MatSolveTrans - Solves A' x = b, given a factored matrix.
1457 
1458    Input Parameters:
1459 .  mat - the factored matrix
1460 .  b - the right-hand-side vector
1461 
1462    Output Parameter:
1463 .  x - the result vector
1464 
1465    Collective on Mat and Vec
1466 
1467    Notes:
1468    The vectors b and x cannot be the same.  I.e., one cannot
1469    call MatSolveTrans(A,x,x).
1470 
1471 .keywords: matrix, linear system, solve, LU, Cholesky, transpose
1472 
1473 .seealso: MatSolve(), MatSolveAdd(), MatSolveTransAdd()
1474 @*/
1475 int MatSolveTrans(Mat mat,Vec b,Vec x)
1476 {
1477   int ierr;
1478 
1479   PetscFunctionBegin;
1480   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1481   PetscValidHeaderSpecific(b,VEC_COOKIE);  PetscValidHeaderSpecific(x,VEC_COOKIE);
1482   if (!mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Unfactored matrix");
1483   if (x == b) SETERRQ(PETSC_ERR_ARG_IDN,0,"x and b must be different vectors");
1484   if (!mat->ops->solvetrans) SETERRQ(PETSC_ERR_SUP,0,"");
1485   if (mat->M != x->N) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec x: global dim");
1486   if (mat->N != b->N) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec b: global dim");
1487 
1488   PLogEventBegin(MAT_SolveTrans,mat,b,x,0);
1489   ierr = (*mat->ops->solvetrans)(mat,b,x); CHKERRQ(ierr);
1490   PLogEventEnd(MAT_SolveTrans,mat,b,x,0);
1491   PetscFunctionReturn(0);
1492 }
1493 
1494 #undef __FUNC__
1495 #define __FUNC__ "MatSolveTransAdd"
1496 /*@
1497    MatSolveTransAdd - Computes x = y + inv(trans(A)) b, given a
1498                       factored matrix.
1499 
1500    Input Parameters:
1501 .  mat - the factored matrix
1502 .  b - the right-hand-side vector
1503 .  y - the vector to be added to
1504 
1505    Output Parameter:
1506 .  x - the result vector
1507 
1508    Collective on Mat and Vec
1509 
1510    Notes:
1511    The vectors b and x cannot be the same.  I.e., one cannot
1512    call MatSolveTransAdd(A,x,y,x).
1513 
1514 .keywords: matrix, linear system, solve, LU, Cholesky, transpose, add
1515 
1516 .seealso: MatSolve(), MatSolveAdd(), MatSolveTrans()
1517 @*/
1518 int MatSolveTransAdd(Mat mat,Vec b,Vec y,Vec x)
1519 {
1520   Scalar one = 1.0;
1521   int    ierr;
1522   Vec    tmp;
1523 
1524   PetscFunctionBegin;
1525   PetscValidHeaderSpecific(mat,MAT_COOKIE);PetscValidHeaderSpecific(y,VEC_COOKIE);
1526   PetscValidHeaderSpecific(b,VEC_COOKIE);  PetscValidHeaderSpecific(x,VEC_COOKIE);
1527   if (x == b) SETERRQ(PETSC_ERR_ARG_IDN,0,"x and b must be different vectors");
1528   if (!mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Unfactored matrix");
1529   if (mat->M != x->N) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec x: global dim");
1530   if (mat->N != b->N) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec b: global dim");
1531   if (mat->N != y->N) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec y: global dim");
1532   if (x->n != y->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Vec x,Vec y: local dim");
1533 
1534   PLogEventBegin(MAT_SolveTransAdd,mat,b,x,y);
1535   if (mat->ops->solvetransadd) {
1536     ierr = (*mat->ops->solvetransadd)(mat,b,y,x); CHKERRQ(ierr);
1537   }
1538   else {
1539     /* do the solve then the add manually */
1540     if (x != y) {
1541       ierr = MatSolveTrans(mat,b,x); CHKERRQ(ierr);
1542       ierr = VecAXPY(&one,y,x); CHKERRQ(ierr);
1543     }
1544     else {
1545       ierr = VecDuplicate(x,&tmp); CHKERRQ(ierr);
1546       PLogObjectParent(mat,tmp);
1547       ierr = VecCopy(x,tmp); CHKERRQ(ierr);
1548       ierr = MatSolveTrans(mat,b,x); CHKERRQ(ierr);
1549       ierr = VecAXPY(&one,tmp,x); CHKERRQ(ierr);
1550       ierr = VecDestroy(tmp); CHKERRQ(ierr);
1551     }
1552   }
1553   PLogEventEnd(MAT_SolveTransAdd,mat,b,x,y);
1554   PetscFunctionReturn(0);
1555 }
1556 /* ----------------------------------------------------------------*/
1557 
1558 #undef __FUNC__
1559 #define __FUNC__ "MatRelax"
1560 /*@
1561    MatRelax - Computes one relaxation sweep.
1562 
1563    Input Parameters:
1564 .  mat - the matrix
1565 .  b - the right hand side
1566 .  omega - the relaxation factor
1567 .  flag - flag indicating the type of SOR, one of
1568 $     SOR_FORWARD_SWEEP
1569 $     SOR_BACKWARD_SWEEP
1570 $     SOR_SYMMETRIC_SWEEP (SSOR method)
1571 $     SOR_LOCAL_FORWARD_SWEEP
1572 $     SOR_LOCAL_BACKWARD_SWEEP
1573 $     SOR_LOCAL_SYMMETRIC_SWEEP (local SSOR)
1574 $     SOR_APPLY_UPPER, SOR_APPLY_LOWER - applies
1575 $       upper/lower triangular part of matrix to
1576 $       vector (with omega)
1577 $     SOR_ZERO_INITIAL_GUESS - zero initial guess
1578 .  shift -  diagonal shift
1579 .  its - the number of iterations
1580 
1581    Output Parameters:
1582 .  x - the solution (can contain an initial guess)
1583 
1584    Collective on Mat and Vec
1585 
1586    Notes:
1587    SOR_LOCAL_FORWARD_SWEEP, SOR_LOCAL_BACKWARD_SWEEP, and
1588    SOR_LOCAL_SYMMETRIC_SWEEP perform seperate independent smoothings
1589    on each processor.
1590 
1591    Application programmers will not generally use MatRelax() directly,
1592    but instead will employ the SLES/PC interface.
1593 
1594    Notes for Advanced Users:
1595    The flags are implemented as bitwise inclusive or operations.
1596    For example, use (SOR_ZERO_INITIAL_GUESS | SOR_SYMMETRIC_SWEEP)
1597    to specify a zero initial guess for SSOR.
1598 
1599 .keywords: matrix, relax, relaxation, sweep
1600 @*/
1601 int MatRelax(Mat mat,Vec b,double omega,MatSORType flag,double shift,
1602              int its,Vec x)
1603 {
1604   int ierr;
1605 
1606   PetscFunctionBegin;
1607   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1608   PetscValidHeaderSpecific(b,VEC_COOKIE);  PetscValidHeaderSpecific(x,VEC_COOKIE);
1609   if (!mat->ops->relax) SETERRQ(PETSC_ERR_SUP,0,"");
1610   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix");
1611   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix");
1612   if (mat->N != x->N) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec x: global dim");
1613   if (mat->M != b->N) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec b: global dim");
1614   if (mat->m != b->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec b: local dim");
1615 
1616   PLogEventBegin(MAT_Relax,mat,b,x,0);
1617   ierr =(*mat->ops->relax)(mat,b,omega,flag,shift,its,x); CHKERRQ(ierr);
1618   PLogEventEnd(MAT_Relax,mat,b,x,0);
1619   PetscFunctionReturn(0);
1620 }
1621 
1622 #undef __FUNC__
1623 #define __FUNC__ "MatCopy_Basic"
1624 /*
1625       Default matrix copy routine.
1626 */
1627 int MatCopy_Basic(Mat A,Mat B)
1628 {
1629   int    ierr,i,rstart,rend,nz,*cwork;
1630   Scalar *vwork;
1631 
1632   PetscFunctionBegin;
1633   ierr = MatZeroEntries(B); CHKERRQ(ierr);
1634   ierr = MatGetOwnershipRange(A,&rstart,&rend); CHKERRQ(ierr);
1635   for (i=rstart; i<rend; i++) {
1636     ierr = MatGetRow(A,i,&nz,&cwork,&vwork); CHKERRQ(ierr);
1637     ierr = MatSetValues(B,1,&i,nz,cwork,vwork,INSERT_VALUES); CHKERRQ(ierr);
1638     ierr = MatRestoreRow(A,i,&nz,&cwork,&vwork); CHKERRQ(ierr);
1639   }
1640   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
1641   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
1642   PetscFunctionReturn(0);
1643 }
1644 
1645 #undef __FUNC__
1646 #define __FUNC__ "MatCopy"
1647 /*@C
1648    MatCopy - Copys a matrix to another matrix.
1649 
1650    Input Parameters:
1651 .  A - the matrix
1652 
1653    Output Parameter:
1654 .  B - where the copy is put
1655 
1656    Collective on Mat
1657 
1658    Notes:
1659    MatCopy() copies the matrix entries of a matrix to another existing
1660    matrix (after first zeroing the second matrix).  A related routine is
1661    MatConvert(), which first creates a new matrix and then copies the data.
1662 
1663 .keywords: matrix, copy, convert
1664 
1665 .seealso: MatConvert()
1666 @*/
1667 int MatCopy(Mat A,Mat B)
1668 {
1669   int ierr;
1670 
1671   PetscFunctionBegin;
1672   PetscValidHeaderSpecific(A,MAT_COOKIE); PetscValidHeaderSpecific(B,MAT_COOKIE);
1673   if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix");
1674   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix");
1675   if (A->M != B->M || A->N != B->N) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Mat A,Mat B: global dim");
1676 
1677   PLogEventBegin(MAT_Copy,A,B,0,0);
1678   if (A->ops->copy) {
1679     ierr = (*A->ops->copy)(A,B); CHKERRQ(ierr);
1680   }
1681   else { /* generic conversion */
1682     ierr = MatCopy_Basic(A,B); CHKERRQ(ierr);
1683   }
1684   PLogEventEnd(MAT_Copy,A,B,0,0);
1685   PetscFunctionReturn(0);
1686 }
1687 
1688 static int MatConvertersSet = 0;
1689 static int (*MatConverters[MAX_MATRIX_TYPES][MAX_MATRIX_TYPES])(Mat,MatType,Mat*) =
1690            {{0,0,0,0,0,0,0,0,0,0,0,0},{0,0,0,0,0,0,0,0,0,0,0,0},
1691             {0,0,0,0,0,0,0,0,0,0,0,0},{0,0,0,0,0,0,0,0,0,0,0,0},
1692             {0,0,0,0,0,0,0,0,0,0,0,0},{0,0,0,0,0,0,0,0,0,0,0,0},
1693             {0,0,0,0,0,0,0,0,0,0,0,0},{0,0,0,0,0,0,0,0,0,0,0,0},
1694             {0,0,0,0,0,0,0,0,0,0,0,0},{0,0,0,0,0,0,0,0,0,0,0,0},
1695             {0,0,0,0,0,0,0,0,0,0,0,0},{0,0,0,0,0,0,0,0,0,0,0,0}};
1696 
1697 #undef __FUNC__
1698 #define __FUNC__ "MatConvertRegister"
1699 /*@C
1700     MatConvertRegister - Allows one to register a routine that converts between
1701         two matrix types.
1702 
1703   Input Parameters:
1704 .   intype - the type of matrix (defined in include/mat.h), for example, MATSEQAIJ.
1705 .   outtype - new matrix type, or MATSAME
1706 
1707   Not Collective
1708 
1709 .seealso: MatConvertRegisterAll()
1710 @*/
1711 int MatConvertRegister(MatType intype,MatType outtype,int (*converter)(Mat,MatType,Mat*))
1712 {
1713   PetscFunctionBegin;
1714   MatConverters[intype][outtype] = converter;
1715   MatConvertersSet               = 1;
1716   PetscFunctionReturn(0);
1717 }
1718 
1719 #undef __FUNC__
1720 #define __FUNC__ "MatConvert"
1721 /*@C
1722    MatConvert - Converts a matrix to another matrix, either of the same
1723    or different type.
1724 
1725    Input Parameters:
1726 .  mat - the matrix
1727 .  newtype - new matrix type.  Use MATSAME to create a new matrix of the
1728    same type as the original matrix.
1729 
1730    Output Parameter:
1731 .  M - pointer to place new matrix
1732 
1733    Collective on Mat
1734 
1735    Notes:
1736    MatConvert() first creates a new matrix and then copies the data from
1737    the first matrix.  A related routine is MatCopy(), which copies the matrix
1738    entries of one matrix to another already existing matrix context.
1739 
1740 .keywords: matrix, copy, convert
1741 
1742 .seealso: MatCopy(), MatDuplicate()
1743 @*/
1744 int MatConvert(Mat mat,MatType newtype,Mat *M)
1745 {
1746   int ierr;
1747 
1748   PetscFunctionBegin;
1749   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1750   PetscValidPointer(M);
1751   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix");
1752   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix");
1753 
1754   if (newtype > MAX_MATRIX_TYPES || newtype < -1) {
1755     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,1,"Not a valid matrix type");
1756   }
1757   *M  = 0;
1758 
1759   if (!MatConvertersSet) {
1760     ierr = MatLoadRegisterAll(); CHKERRQ(ierr);
1761   }
1762 
1763   PLogEventBegin(MAT_Convert,mat,0,0,0);
1764   if ((newtype == mat->type || newtype == MATSAME) && mat->ops->convertsametype) {
1765     ierr = (*mat->ops->convertsametype)(mat,M,COPY_VALUES); CHKERRQ(ierr);
1766   } else {
1767     if (!MatConvertersSet) {
1768       ierr = MatConvertRegisterAll(); CHKERRQ(ierr);
1769     }
1770     if (!MatConverters[mat->type][newtype]) {
1771       SETERRQ(PETSC_ERR_ARG_WRONG,1,"Invalid matrix type, or matrix converter not registered");
1772     }
1773     ierr = (*MatConverters[mat->type][newtype])(mat,newtype,M); CHKERRQ(ierr);
1774   }
1775   PLogEventEnd(MAT_Convert,mat,0,0,0);
1776   PetscFunctionReturn(0);
1777 }
1778 
1779 #undef __FUNC__
1780 #define __FUNC__ "MatDuplicate"
1781 /*@C
1782    MatDuplicate - Duplicates a matrix including the non-zero structure, but
1783      does not copy over the values.
1784 
1785    Input Parameters:
1786 .  mat - the matrix
1787 
1788    Output Parameter:
1789 .  M - pointer to place new matrix
1790 
1791    Collective on Mat
1792 
1793 .keywords: matrix, copy, convert, duplicate
1794 
1795 .seealso: MatCopy(), MatDuplicate(), MatConvert()
1796 @*/
1797 int MatDuplicate(Mat mat,Mat *M)
1798 {
1799   int ierr;
1800 
1801   PetscFunctionBegin;
1802   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1803   PetscValidPointer(M);
1804   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix");
1805   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix");
1806 
1807   *M  = 0;
1808   PLogEventBegin(MAT_Convert,mat,0,0,0);
1809   if (!mat->ops->convertsametype) {
1810     SETERRQ(PETSC_ERR_SUP,1,"Not written for this matrix type");
1811   }
1812   ierr = (*mat->ops->convertsametype)(mat,M,DO_NOT_COPY_VALUES); CHKERRQ(ierr);
1813   PLogEventEnd(MAT_Convert,mat,0,0,0);
1814   PetscFunctionReturn(0);
1815 }
1816 
1817 #undef __FUNC__
1818 #define __FUNC__ "MatGetDiagonal"
1819 /*@
1820    MatGetDiagonal - Gets the diagonal of a matrix.
1821 
1822    Input Parameters:
1823 .  mat - the matrix
1824 .  v - the vector for storing the diagonal
1825 
1826    Output Parameter:
1827 .  v - the diagonal of the matrix
1828 
1829    Collective on Mat and Vec
1830 
1831    Notes: For the SeqAIJ matrix format, this routine may also be called
1832     on a LU factored matrix; in that case it routines the reciprocal of
1833     the diagonal entries in U. It returns the entries permuted by the
1834     row and column permutation used during the symbolic factorization.
1835 
1836 .keywords: matrix, get, diagonal
1837 @*/
1838 int MatGetDiagonal(Mat mat,Vec v)
1839 {
1840   int ierr;
1841 
1842   PetscFunctionBegin;
1843   PetscValidHeaderSpecific(mat,MAT_COOKIE);PetscValidHeaderSpecific(v,VEC_COOKIE);
1844   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix");
1845   /*
1846      The error checking for a factored matrix is handled inside
1847     each matrix type, since MatGetDiagonal() is supported by
1848     factored AIJ matrices
1849   */
1850   /* if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix");  */
1851   if (!mat->ops->getdiagonal) SETERRQ(PETSC_ERR_SUP,0,"");
1852   ierr = (*mat->ops->getdiagonal)(mat,v);CHKERRQ(ierr);
1853   PetscFunctionReturn(0);
1854 }
1855 
1856 #undef __FUNC__
1857 #define __FUNC__ "MatTranspose"
1858 /*@C
1859    MatTranspose - Computes an in-place or out-of-place transpose of a matrix.
1860 
1861    Input Parameter:
1862 .  mat - the matrix to transpose
1863 
1864    Output Parameters:
1865 .  B - the transpose (or pass in PETSC_NULL for an in-place transpose)
1866 
1867    Collective on Mat
1868 
1869 .keywords: matrix, transpose
1870 
1871 .seealso: MatMultTrans(), MatMultTransAdd()
1872 @*/
1873 int MatTranspose(Mat mat,Mat *B)
1874 {
1875   int ierr;
1876 
1877   PetscFunctionBegin;
1878   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1879   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix");
1880   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix");
1881   if (!mat->ops->transpose) SETERRQ(PETSC_ERR_SUP,0,"");
1882   ierr = (*mat->ops->transpose)(mat,B);CHKERRQ(ierr);
1883   PetscFunctionReturn(0);
1884 }
1885 
1886 #undef __FUNC__
1887 #define __FUNC__ "MatPermute"
1888 /*@C
1889    MatPermute - Creates a new matrix with rows and columns permuted from the
1890        original.
1891 
1892    Input Parameter:
1893 .  mat - the matrix to permute
1894 .  row - row permutation
1895 .  col - column permutation
1896 
1897    Output Parameters:
1898 .  B - the permuted matrix
1899 
1900    Collective on Mat
1901 
1902 .keywords: matrix, transpose
1903 
1904 .seealso: MatGetReordering()
1905 @*/
1906 int MatPermute(Mat mat,IS row,IS col,Mat *B)
1907 {
1908   int ierr;
1909 
1910   PetscFunctionBegin;
1911   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1912   PetscValidHeaderSpecific(row,IS_COOKIE);
1913   PetscValidHeaderSpecific(col,IS_COOKIE);
1914   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix");
1915   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix");
1916   if (!mat->ops->permute) SETERRQ(PETSC_ERR_SUP,0,"");
1917   ierr = (*mat->ops->permute)(mat,row,col,B);CHKERRQ(ierr);
1918   PetscFunctionReturn(0);
1919 }
1920 
1921 #undef __FUNC__
1922 #define __FUNC__ "MatEqual"
1923 /*@
1924    MatEqual - Compares two matrices.
1925 
1926    Input Parameters:
1927 .  A - the first matrix
1928 .  B - the second matrix
1929 
1930    Output Parameter:
1931 .  flg : PETSC_TRUE if the matrices are equal;
1932          PETSC_FALSE otherwise.
1933 
1934    Collective on Mat
1935 
1936 .keywords: matrix, equal, equivalent
1937 @*/
1938 int MatEqual(Mat A,Mat B,PetscTruth *flg)
1939 {
1940   int ierr;
1941 
1942   PetscFunctionBegin;
1943   PetscValidHeaderSpecific(A,MAT_COOKIE); PetscValidHeaderSpecific(B,MAT_COOKIE);
1944   PetscValidIntPointer(flg);
1945   if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix");
1946   if (!B->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix");
1947   if (A->M != B->M || A->N != B->N) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Mat A,Mat B: global dim");
1948   if (!A->ops->equal) SETERRQ(PETSC_ERR_SUP,0,"");
1949   ierr = (*A->ops->equal)(A,B,flg);CHKERRQ(ierr);
1950   PetscFunctionReturn(0);
1951 }
1952 
1953 #undef __FUNC__
1954 #define __FUNC__ "MatDiagonalScale"
1955 /*@
1956    MatDiagonalScale - Scales a matrix on the left and right by diagonal
1957    matrices that are stored as vectors.  Either of the two scaling
1958    matrices can be PETSC_NULL.
1959 
1960    Input Parameters:
1961 .  mat - the matrix to be scaled
1962 .  l - the left scaling vector (or PETSC_NULL)
1963 .  r - the right scaling vector (or PETSC_NULL)
1964 
1965    Collective on Mat
1966 
1967    Notes:
1968    MatDiagonalScale() computes A <- LAR, where
1969 $      L = a diagonal matrix
1970 $      R = a diagonal matrix
1971 
1972 .keywords: matrix, diagonal, scale
1973 
1974 .seealso: MatDiagonalScale()
1975 @*/
1976 int MatDiagonalScale(Mat mat,Vec l,Vec r)
1977 {
1978   int ierr;
1979 
1980   PetscFunctionBegin;
1981   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1982   if (!mat->ops->diagonalscale) SETERRQ(PETSC_ERR_SUP,0,"");
1983   if (l) PetscValidHeaderSpecific(l,VEC_COOKIE);
1984   if (r) PetscValidHeaderSpecific(r,VEC_COOKIE);
1985   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix");
1986   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix");
1987 
1988   PLogEventBegin(MAT_Scale,mat,0,0,0);
1989   ierr = (*mat->ops->diagonalscale)(mat,l,r); CHKERRQ(ierr);
1990   PLogEventEnd(MAT_Scale,mat,0,0,0);
1991   PetscFunctionReturn(0);
1992 }
1993 
1994 #undef __FUNC__
1995 #define __FUNC__ "MatScale"
1996 /*@
1997     MatScale - Scales all elements of a matrix by a given number.
1998 
1999     Input Parameters:
2000 .   mat - the matrix to be scaled
2001 .   a  - the scaling value
2002 
2003     Output Parameter:
2004 .   mat - the scaled matrix
2005 
2006    Collective on Mat
2007 
2008 .keywords: matrix, scale
2009 
2010 .seealso: MatDiagonalScale()
2011 @*/
2012 int MatScale(Scalar *a,Mat mat)
2013 {
2014   int ierr;
2015 
2016   PetscFunctionBegin;
2017   PetscValidHeaderSpecific(mat,MAT_COOKIE);
2018   PetscValidScalarPointer(a);
2019   if (!mat->ops->scale) SETERRQ(PETSC_ERR_SUP,0,"");
2020   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix");
2021   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix");
2022 
2023   PLogEventBegin(MAT_Scale,mat,0,0,0);
2024   ierr = (*mat->ops->scale)(a,mat); CHKERRQ(ierr);
2025   PLogEventEnd(MAT_Scale,mat,0,0,0);
2026   PetscFunctionReturn(0);
2027 }
2028 
2029 #undef __FUNC__
2030 #define __FUNC__ "MatNorm"
2031 /*@
2032    MatNorm - Calculates various norms of a matrix.
2033 
2034    Input Parameters:
2035 .  mat - the matrix
2036 .  type - the type of norm, NORM_1, NORM_2, NORM_FROBENIUS, NORM_INFINITY
2037 
2038    Output Parameters:
2039 .  norm - the resulting norm
2040 
2041    Collective on Mat
2042 
2043 .keywords: matrix, norm, Frobenius
2044 @*/
2045 int MatNorm(Mat mat,NormType type,double *norm)
2046 {
2047   int ierr;
2048 
2049   PetscFunctionBegin;
2050   PetscValidHeaderSpecific(mat,MAT_COOKIE);
2051   PetscValidScalarPointer(norm);
2052 
2053   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix");
2054   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix");
2055   if (!mat->ops->norm) SETERRQ(PETSC_ERR_SUP,0,"Not for this matrix type");
2056   ierr = (*mat->ops->norm)(mat,type,norm);CHKERRQ(ierr);
2057   PetscFunctionReturn(0);
2058 }
2059 
2060 /*
2061      This variable is used to prevent counting of MatAssemblyBegin() that
2062    are called from within a MatAssemblyEnd().
2063 */
2064 static int MatAssemblyEnd_InUse = 0;
2065 #undef __FUNC__
2066 #define __FUNC__ "MatAssemblyBegin"
2067 /*@
2068    MatAssemblyBegin - Begins assembling the matrix.  This routine should
2069    be called after completing all calls to MatSetValues().
2070 
2071    Input Parameters:
2072 .  mat - the matrix
2073 .  type - type of assembly, either MAT_FLUSH_ASSEMBLY or MAT_FINAL_ASSEMBLY
2074 
2075    Collective on Mat
2076 
2077    Notes:
2078    MatSetValues() generally caches the values.  The matrix is ready to
2079    use only after MatAssemblyBegin() and MatAssemblyEnd() have been called.
2080    Use MAT_FLUSH_ASSEMBLY when switching between ADD_VALUES and INSERT_VALUES
2081    in MatSetValues(); use MAT_FINAL_ASSEMBLY for the final assembly before
2082    using the matrix.
2083 
2084 .keywords: matrix, assembly, assemble, begin
2085 
2086 .seealso: MatAssemblyEnd(), MatSetValues()
2087 @*/
2088 int MatAssemblyBegin(Mat mat,MatAssemblyType type)
2089 {
2090   int ierr;
2091 
2092   PetscFunctionBegin;
2093   PetscValidHeaderSpecific(mat,MAT_COOKIE);
2094   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix.\n did you forget to call MatSetUnfactored()?");
2095   if (mat->assembled) {
2096     mat->was_assembled = PETSC_TRUE;
2097     mat->assembled     = PETSC_FALSE;
2098   }
2099   if (!MatAssemblyEnd_InUse) {
2100     PLogEventBegin(MAT_AssemblyBegin,mat,0,0,0);
2101     if (mat->ops->assemblybegin){ierr = (*mat->ops->assemblybegin)(mat,type);CHKERRQ(ierr);}
2102     PLogEventEnd(MAT_AssemblyBegin,mat,0,0,0);
2103   } else {
2104     if (mat->ops->assemblybegin){ierr = (*mat->ops->assemblybegin)(mat,type);CHKERRQ(ierr);}
2105   }
2106   PetscFunctionReturn(0);
2107 }
2108 
2109 
2110 #undef __FUNC__
2111 #define __FUNC__ "MatView_Private"
2112 /*
2113     Processes command line options to determine if/how a matrix
2114   is to be viewed. Called by MatAssemblyEnd() and MatLoad().
2115 */
2116 int MatView_Private(Mat mat)
2117 {
2118   int ierr,flg;
2119 
2120   PetscFunctionBegin;
2121   ierr = OptionsHasName(PETSC_NULL,"-mat_view_info",&flg); CHKERRQ(ierr);
2122   if (flg) {
2123     ierr = ViewerPushFormat(VIEWER_STDOUT_(mat->comm),VIEWER_FORMAT_ASCII_INFO,0);CHKERRQ(ierr);
2124     ierr = MatView(mat,VIEWER_STDOUT_(mat->comm)); CHKERRQ(ierr);
2125     ierr = ViewerPopFormat(VIEWER_STDOUT_(mat->comm));CHKERRQ(ierr);
2126   }
2127   ierr = OptionsHasName(PETSC_NULL,"-mat_view_info_detailed",&flg);CHKERRQ(ierr);
2128   if (flg) {
2129     ierr = ViewerPushFormat(VIEWER_STDOUT_(mat->comm),VIEWER_FORMAT_ASCII_INFO_LONG,0);CHKERRQ(ierr);
2130     ierr = MatView(mat,VIEWER_STDOUT_(mat->comm)); CHKERRQ(ierr);
2131     ierr = ViewerPopFormat(VIEWER_STDOUT_(mat->comm));CHKERRQ(ierr);
2132   }
2133   ierr = OptionsHasName(PETSC_NULL,"-mat_view",&flg); CHKERRQ(ierr);
2134   if (flg) {
2135     ierr = MatView(mat,VIEWER_STDOUT_(mat->comm)); CHKERRQ(ierr);
2136   }
2137   ierr = OptionsHasName(PETSC_NULL,"-mat_view_matlab",&flg); CHKERRQ(ierr);
2138   if (flg) {
2139     ierr = ViewerPushFormat(VIEWER_STDOUT_(mat->comm),VIEWER_FORMAT_ASCII_MATLAB,"M");CHKERRQ(ierr);
2140     ierr = MatView(mat,VIEWER_STDOUT_(mat->comm)); CHKERRQ(ierr);
2141     ierr = ViewerPopFormat(VIEWER_STDOUT_(mat->comm));CHKERRQ(ierr);
2142   }
2143   ierr = OptionsHasName(PETSC_NULL,"-mat_view_draw",&flg); CHKERRQ(ierr);
2144   if (flg) {
2145     ierr = OptionsHasName(0,"-mat_view_contour",&flg); CHKERRQ(ierr);
2146     if (flg) {
2147       ViewerPushFormat(VIEWER_DRAWX_(mat->comm),VIEWER_FORMAT_DRAW_CONTOUR,0);CHKERRQ(ierr);
2148     }
2149     ierr = MatView(mat,VIEWER_DRAWX_(mat->comm)); CHKERRQ(ierr);
2150     ierr = ViewerFlush(VIEWER_DRAWX_(mat->comm)); CHKERRQ(ierr);
2151     if (flg) {
2152       ViewerPopFormat(VIEWER_DRAWX_(mat->comm));CHKERRQ(ierr);
2153     }
2154   }
2155   PetscFunctionReturn(0);
2156 }
2157 
2158 #undef __FUNC__
2159 #define __FUNC__ "MatAssemblyEnd"
2160 /*@
2161    MatAssemblyEnd - Completes assembling the matrix.  This routine should
2162    be called after MatAssemblyBegin().
2163 
2164    Input Parameters:
2165 .  mat - the matrix
2166 .  type - type of assembly, either MAT_FLUSH_ASSEMBLY or MAT_FINAL_ASSEMBLY
2167 
2168    Collective on Mat
2169 
2170    Options Database Keys:
2171 $  -mat_view_info : Prints info on matrix at
2172 $      conclusion of MatEndAssembly()
2173 $  -mat_view_info_detailed: Prints more detailed info.
2174 $  -mat_view : Prints matrix in ASCII format.
2175 $  -mat_view_matlab : Prints matrix in Matlab format.
2176 $  -mat_view_draw : Draws nonzero structure of matrix,
2177 $      using MatView() and DrawOpenX().
2178 $  -display <name> : Set display name (default is host)
2179 $  -draw_pause <sec> : Set number of seconds to pause after display
2180 
2181    Notes:
2182    MatSetValues() generally caches the values.  The matrix is ready to
2183    use only after MatAssemblyBegin() and MatAssemblyEnd() have been called.
2184    Use MAT_FLUSH_ASSEMBLY when switching between ADD_VALUES and INSERT_VALUES
2185    in MatSetValues(); use MAT_FINAL_ASSEMBLY for the final assembly before
2186    using the matrix.
2187 
2188 .keywords: matrix, assembly, assemble, end
2189 
2190 .seealso: MatAssemblyBegin(), MatSetValues(), DrawOpenX(), MatView()
2191 @*/
2192 int MatAssemblyEnd(Mat mat,MatAssemblyType type)
2193 {
2194   int        ierr;
2195   static int inassm = 0;
2196 
2197   PetscFunctionBegin;
2198   PetscValidHeaderSpecific(mat,MAT_COOKIE);
2199 
2200   inassm++;
2201   MatAssemblyEnd_InUse++;
2202   PLogEventBegin(MAT_AssemblyEnd,mat,0,0,0);
2203   if (mat->ops->assemblyend) {
2204     ierr = (*mat->ops->assemblyend)(mat,type); CHKERRQ(ierr);
2205   }
2206 
2207   /* Flush assembly is not a true assembly */
2208   if (type != MAT_FLUSH_ASSEMBLY) {
2209     mat->assembled  = PETSC_TRUE; mat->num_ass++;
2210   }
2211   mat->insertmode = NOT_SET_VALUES;
2212   PLogEventEnd(MAT_AssemblyEnd,mat,0,0,0);
2213   MatAssemblyEnd_InUse--;
2214 
2215   if (inassm == 1 && type != MAT_FLUSH_ASSEMBLY) {
2216     ierr = MatView_Private(mat); CHKERRQ(ierr);
2217   }
2218   inassm--;
2219   PetscFunctionReturn(0);
2220 }
2221 
2222 
2223 #undef __FUNC__
2224 #define __FUNC__ "MatCompress"
2225 /*@
2226    MatCompress - Tries to store the matrix in as little space as
2227    possible.  May fail if memory is already fully used, since it
2228    tries to allocate new space.
2229 
2230    Input Parameters:
2231 .  mat - the matrix
2232 
2233    Collective on Mat
2234 
2235 .keywords: matrix, compress
2236 @*/
2237 int MatCompress(Mat mat)
2238 {
2239   int ierr;
2240 
2241   PetscFunctionBegin;
2242   PetscValidHeaderSpecific(mat,MAT_COOKIE);
2243   if (mat->ops->compress) {ierr = (*mat->ops->compress)(mat);CHKERRQ(ierr);}
2244   PetscFunctionReturn(0);
2245 }
2246 
2247 #undef __FUNC__
2248 #define __FUNC__ "MatSetOption"
2249 /*@
2250    MatSetOption - Sets a parameter option for a matrix. Some options
2251    may be specific to certain storage formats.  Some options
2252    determine how values will be inserted (or added). Sorted,
2253    row-oriented input will generally assemble the fastest. The default
2254    is row-oriented, nonsorted input.
2255 
2256    Input Parameters:
2257 .  mat - the matrix
2258 .  option - the option, one of those listed below (and possibly others),
2259              e.g., MAT_ROWS_SORTED, MAT_NEW_NONZERO_LOCATION_ERROR
2260 
2261    Collective on Mat
2262 
2263    Options Describing Matrix Structure:
2264 $    MAT_SYMMETRIC - symmetric in terms of both structure and value
2265 $    MAT_STRUCTURALLY_SYMMETRIC - symmetric nonzero structure
2266 
2267    Options For Use with MatSetValues():
2268    Insert a logically dense subblock, which can be
2269 
2270 $    MAT_ROW_ORIENTED - row-oriented
2271 $    MAT_COLUMN_ORIENTED - column-oriented
2272 $    MAT_ROWS_SORTED - sorted by row
2273 $    MAT_ROWS_UNSORTED - not sorted by row
2274 $    MAT_COLUMNS_SORTED - sorted by column
2275 $    MAT_COLUMNS_UNSORTED - not sorted by column
2276 
2277    Not these options reflect the data you pass in with MatSetValues(); it has
2278    nothing to do with how the data is stored internally in the matrix
2279    data structure.
2280 
2281    When (re)assembling a matrix, we can restrict the input for
2282    efficiency/debugging purposes.
2283 
2284 $    MAT_NO_NEW_NONZERO_LOCATIONS - additional insertions will not be
2285         allowed if they generate a new nonzero
2286 $    MAT_YES_NEW_NONZERO_LOCATIONS - additional insertions will be allowed
2287 $    MAT_NO_NEW_DIAGONALS - additional insertions will not be allowed if
2288          they generate a nonzero in a new diagonal (for block diagonal format only)
2289 $    MAT_YES_NEW_DIAGONALS - new diagonals will be allowed (for block diagonal format only)
2290 $    MAT_IGNORE_OFF_PROC_ENTRIES - drop off-processor entries
2291 $    MAT_NEW_NONZERO_LOCATION_ERROR - generate error for new matrix entry
2292 $    MAT_USE_HASH_TABLE - use hash table which speeds up the Matrix assembly
2293 
2294    Notes:
2295    Some options are relevant only for particular matrix types and
2296    are thus ignored by others.  Other options are not supported by
2297    certain matrix types and will generate an error message if set.
2298 
2299    If using a Fortran 77 module to compute a matrix, one may need to
2300    use the column-oriented option (or convert to the row-oriented
2301    format).
2302 
2303    MAT_NO_NEW_NONZERO_LOCATIONS indicates that any add or insertion
2304    that would generate a new entry in the nonzero structure is instead
2305    ignored.  Thus, if memory has not alredy been allocated for this particular
2306    data, then the insertion is ignored. For dense matrices, in which
2307    the entire array is allocated, no entries are ever ignored.
2308 
2309    MAT_NEW_NONZERO_LOCATION_ERROR indicates that any add or insertion
2310    that would generate a new entry in the nonzero structure instead produces
2311    an error. (Currently supported for AIJ and BAIJ formats only.)
2312    This is a useful flag when using SAME_NONZERO_PATTERN in calling
2313    SLESSetOperators() to ensure that the nonzero pattern truely does
2314    remain unchanged.
2315 
2316    MAT_NEW_NONZERO_ALLOCATION_ERROR indicates that any add or insertion
2317    that would generate a new entry that has not been preallocated will
2318    instead produce an error. (Currently supported for AIJ and BAIJ formats
2319    only.) This is a useful flag when debugging matrix memory preallocation.
2320 
2321    MAT_IGNORE_OFF_PROC_ENTRIES indicates entries destined for
2322    other processors should be dropped, rather than stashed.
2323    This is useful if you know that the "owning" processor is also
2324    always generating the correct matrix entries, so that PETSc need
2325    not transfer duplicate entries generated on another processor.
2326 
2327    MAT_USE_HASH_TABLE indicates that a hash table be used to improve the
2328    searches during matrix assembly. When this flag is set, the hash table
2329    is created during the first Matrix Assembly. This hash table is
2330    used the next time through, during MatSetVaules()/MatSetVaulesBlocked()
2331    to improve the searching of indices. MAT_NO_NEW_NONZERO_LOCATIONS flag
2332    should be used with MAT_USE_HASH_TABLE flag. This option is currently
2333    supported by MATMPIBAIJ format only.
2334 
2335 .keywords: matrix, option, row-oriented, column-oriented, sorted, nonzero
2336 @*/
2337 int MatSetOption(Mat mat,MatOption op)
2338 {
2339   int ierr;
2340 
2341   PetscFunctionBegin;
2342   PetscValidHeaderSpecific(mat,MAT_COOKIE);
2343   if (mat->ops->setoption) {ierr = (*mat->ops->setoption)(mat,op);CHKERRQ(ierr);}
2344   PetscFunctionReturn(0);
2345 }
2346 
2347 #undef __FUNC__
2348 #define __FUNC__ "MatZeroEntries"
2349 /*@
2350    MatZeroEntries - Zeros all entries of a matrix.  For sparse matrices
2351    this routine retains the old nonzero structure.
2352 
2353    Input Parameters:
2354 .  mat - the matrix
2355 
2356    Collective on Mat
2357 
2358 .keywords: matrix, zero, entries
2359 
2360 .seealso: MatZeroRows()
2361 @*/
2362 int MatZeroEntries(Mat mat)
2363 {
2364   int ierr;
2365 
2366   PetscFunctionBegin;
2367   PetscValidHeaderSpecific(mat,MAT_COOKIE);
2368   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix");
2369   if (!mat->ops->zeroentries) SETERRQ(PETSC_ERR_SUP,0,"");
2370 
2371   PLogEventBegin(MAT_ZeroEntries,mat,0,0,0);
2372   ierr = (*mat->ops->zeroentries)(mat); CHKERRQ(ierr);
2373   PLogEventEnd(MAT_ZeroEntries,mat,0,0,0);
2374   PetscFunctionReturn(0);
2375 }
2376 
2377 #undef __FUNC__
2378 #define __FUNC__ "MatZeroRows"
2379 /*@
2380    MatZeroRows - Zeros all entries (except possibly the main diagonal)
2381    of a set of rows of a matrix.
2382 
2383    Input Parameters:
2384 .  mat - the matrix
2385 .  is - index set of rows to remove
2386 .  diag - pointer to value put in all diagonals of eliminated rows.
2387           Note that diag is not a pointer to an array, but merely a
2388           pointer to a single value.
2389 
2390    Collective on Mat
2391 
2392    Notes:
2393    For the AIJ matrix formats this removes the old nonzero structure,
2394    but does not release memory.  For the dense and block diagonal
2395    formats this does not alter the nonzero structure.
2396 
2397    The user can set a value in the diagonal entry (or for the AIJ and
2398    row formats can optionally remove the main diagonal entry from the
2399    nonzero structure as well, by passing a null pointer as the final
2400    argument).
2401 
2402    For the parallel case, all processes that share the matrix (i.e.,
2403    those in the communicator used for matrix creation) MUST call this
2404    routine, regardless of whether any rows being zeroed are owned by
2405    them.
2406 
2407 .keywords: matrix, zero, rows, boundary conditions
2408 
2409 .seealso: MatZeroEntries(),
2410 @*/
2411 int MatZeroRows(Mat mat,IS is, Scalar *diag)
2412 {
2413   int ierr;
2414 
2415   PetscFunctionBegin;
2416   PetscValidHeaderSpecific(mat,MAT_COOKIE);
2417   PetscValidHeaderSpecific(is,IS_COOKIE);
2418   if (diag) PetscValidScalarPointer(diag);
2419   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix");
2420   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix");
2421   if (!mat->ops->zerorows) SETERRQ(PETSC_ERR_SUP,0,"");
2422 
2423   ierr = (*mat->ops->zerorows)(mat,is,diag); CHKERRQ(ierr);
2424   ierr = MatView_Private(mat); CHKERRQ(ierr);
2425   PetscFunctionReturn(0);
2426 }
2427 
2428 #undef __FUNC__
2429 #define __FUNC__ "MatZeroRowsLocal"
2430 /*@
2431    MatZeroRowsLocal - Zeros all entries (except possibly the main diagonal)
2432    of a set of rows of a matrix; using local numbering of rows.
2433 
2434    Input Parameters:
2435 .  mat - the matrix
2436 .  is - index set of rows to remove
2437 .  diag - pointer to value put in all diagonals of eliminated rows.
2438           Note that diag is not a pointer to an array, but merely a
2439           pointer to a single value.
2440 
2441    Collective on Mat
2442 
2443    Notes:
2444    For the AIJ matrix formats this removes the old nonzero structure,
2445    but does not release memory.  For the dense and block diagonal
2446    formats this does not alter the nonzero structure.
2447 
2448    The user can set a value in the diagonal entry (or for the AIJ and
2449    row formats can optionally remove the main diagonal entry from the
2450    nonzero structure as well, by passing a null pointer as the final
2451    argument).
2452 
2453 .keywords: matrix, zero, rows, boundary conditions
2454 
2455 .seealso: MatZeroEntries(),
2456 @*/
2457 int MatZeroRowsLocal(Mat mat,IS is, Scalar *diag)
2458 {
2459   int ierr;
2460   IS  newis;
2461 
2462   PetscFunctionBegin;
2463   PetscValidHeaderSpecific(mat,MAT_COOKIE);
2464   PetscValidHeaderSpecific(is,IS_COOKIE);
2465   if (diag) PetscValidScalarPointer(diag);
2466   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix");
2467   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix");
2468   if (!mat->ops->zerorows) SETERRQ(PETSC_ERR_SUP,0,"");
2469 
2470   ierr = ISLocalToGlobalMappingApplyIS(mat->mapping,is,&newis); CHKERRQ(ierr);
2471   ierr =  (*mat->ops->zerorows)(mat,newis,diag); CHKERRQ(ierr);
2472   ierr = ISDestroy(newis);
2473   PetscFunctionReturn(0);
2474 }
2475 
2476 #undef __FUNC__
2477 #define __FUNC__ "MatGetSize"
2478 /*@
2479    MatGetSize - Returns the numbers of rows and columns in a matrix.
2480 
2481    Input Parameter:
2482 .  mat - the matrix
2483 
2484    Output Parameters:
2485 .  m - the number of global rows
2486 .  n - the number of global columns
2487 
2488    Not Collective
2489 
2490 .keywords: matrix, dimension, size, rows, columns, global, get
2491 
2492 .seealso: MatGetLocalSize()
2493 @*/
2494 int MatGetSize(Mat mat,int *m,int* n)
2495 {
2496   int ierr;
2497 
2498   PetscFunctionBegin;
2499   PetscValidHeaderSpecific(mat,MAT_COOKIE);
2500   ierr = (*mat->ops->getsize)(mat,m,n);CHKERRQ(ierr);
2501   PetscFunctionReturn(0);
2502 }
2503 
2504 #undef __FUNC__
2505 #define __FUNC__ "MatGetLocalSize"
2506 /*@
2507    MatGetLocalSize - Returns the number of rows and columns in a matrix
2508    stored locally.  This information may be implementation dependent, so
2509    use with care.
2510 
2511    Input Parameters:
2512 .  mat - the matrix
2513 
2514    Output Parameters:
2515 .  m - the number of local rows
2516 .  n - the number of local columns
2517 
2518    Not Collective
2519 
2520 .keywords: matrix, dimension, size, local, rows, columns, get
2521 
2522 .seealso: MatGetSize()
2523 @*/
2524 int MatGetLocalSize(Mat mat,int *m,int* n)
2525 {
2526   int ierr;
2527 
2528   PetscFunctionBegin;
2529   PetscValidHeaderSpecific(mat,MAT_COOKIE);
2530   ierr = (*mat->ops->getlocalsize)(mat,m,n);CHKERRQ(ierr);
2531   PetscFunctionReturn(0);
2532 }
2533 
2534 #undef __FUNC__
2535 #define __FUNC__ "MatGetOwnershipRange"
2536 /*@
2537    MatGetOwnershipRange - Returns the range of matrix rows owned by
2538    this processor, assuming that the matrix is laid out with the first
2539    n1 rows on the first processor, the next n2 rows on the second, etc.
2540    For certain parallel layouts this range may not be well defined.
2541 
2542    Input Parameters:
2543 .  mat - the matrix
2544 
2545    Output Parameters:
2546 .  m - the global index of the first local row
2547 .  n - one more than the global index of the last local row
2548 
2549    Not Collective
2550 
2551 .keywords: matrix, get, range, ownership
2552 @*/
2553 int MatGetOwnershipRange(Mat mat,int *m,int* n)
2554 {
2555   int ierr;
2556 
2557   PetscFunctionBegin;
2558   PetscValidHeaderSpecific(mat,MAT_COOKIE);
2559   PetscValidIntPointer(m);
2560   PetscValidIntPointer(n);
2561   if (!mat->ops->getownershiprange) SETERRQ(PETSC_ERR_SUP,0,"");
2562   ierr = (*mat->ops->getownershiprange)(mat,m,n);CHKERRQ(ierr);
2563   PetscFunctionReturn(0);
2564 }
2565 
2566 #undef __FUNC__
2567 #define __FUNC__ "MatILUFactorSymbolic"
2568 /*@
2569    MatILUFactorSymbolic - Performs symbolic ILU factorization of a matrix.
2570    Uses levels of fill only, not drop tolerance. Use MatLUFactorNumeric()
2571    to complete the factorization.
2572 
2573    Input Parameters:
2574 .  mat - the matrix
2575 .  row - row permutation
2576 .  column - column permutation
2577 .  fill - number of levels of fill
2578 .  f - expected fill as ratio of the original number of nonzeros,
2579        for example 3.0; choosing this parameter well can result in
2580        more efficient use of time and space. Run your code with -log_info
2581        to determine an optimal value to use.
2582 
2583    Output Parameters:
2584 .  fact - new matrix that has been symbolically factored
2585 
2586    Collective on Mat
2587 
2588    Notes:
2589    See the file ${PETSC_DIR}/Performace for additional information about
2590    choosing the fill factor for better efficiency.
2591 
2592 .keywords: matrix, factor, incomplete, ILU, symbolic, fill
2593 
2594 .seealso: MatLUFactorSymbolic(), MatLUFactorNumeric()
2595 @*/
2596 int MatILUFactorSymbolic(Mat mat,IS row,IS col,double f,int fill,Mat *fact)
2597 {
2598   int ierr;
2599 
2600   PetscFunctionBegin;
2601   PetscValidHeaderSpecific(mat,MAT_COOKIE);
2602   PetscValidPointer(fact);
2603   if (fill < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Levels of fill negative");
2604   if (!mat->ops->ilufactorsymbolic) SETERRQ(PETSC_ERR_SUP,0,"Only MatCreateMPIRowbs() matrices support parallel ILU");
2605   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix");
2606   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix");
2607 
2608   PLogEventBegin(MAT_ILUFactorSymbolic,mat,row,col,0);
2609   ierr = (*mat->ops->ilufactorsymbolic)(mat,row,col,f,fill,fact); CHKERRQ(ierr);
2610   PLogEventEnd(MAT_ILUFactorSymbolic,mat,row,col,0);
2611   PetscFunctionReturn(0);
2612 }
2613 
2614 #undef __FUNC__
2615 #define __FUNC__ "MatIncompleteCholeskyFactorSymbolic"
2616 /*@
2617    MatIncompleteCholeskyFactorSymbolic - Performs symbolic incomplete
2618    Cholesky factorization for a symmetric matrix.  Use
2619    MatCholeskyFactorNumeric() to complete the factorization.
2620 
2621    Input Parameters:
2622 .  mat - the matrix
2623 .  perm - row and column permutation
2624 .  fill - levels of fill
2625 .  f - expected fill as ratio of original fill
2626 
2627    Output Parameter:
2628 .  fact - the factored matrix
2629 
2630    Collective on Mat
2631 
2632    Note:  Currently only no-fill factorization is supported.
2633 
2634 .keywords: matrix, factor, incomplete, ICC, Cholesky, symbolic, fill
2635 
2636 .seealso: MatCholeskyFactorNumeric(), MatCholeskyFactor()
2637 @*/
2638 int MatIncompleteCholeskyFactorSymbolic(Mat mat,IS perm,double f,int fill,Mat *fact)
2639 {
2640   int ierr;
2641 
2642   PetscFunctionBegin;
2643   PetscValidHeaderSpecific(mat,MAT_COOKIE);
2644   PetscValidPointer(fact);
2645   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix");
2646   if (fill < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Fill negative");
2647   if (!mat->ops->incompletecholeskyfactorsymbolic) SETERRQ(PETSC_ERR_SUP,0,"Currently only MatCreateMPIRowbs() matrices support ICC in parallel");
2648   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix");
2649 
2650   PLogEventBegin(MAT_IncompleteCholeskyFactorSymbolic,mat,perm,0,0);
2651   ierr = (*mat->ops->incompletecholeskyfactorsymbolic)(mat,perm,f,fill,fact);CHKERRQ(ierr);
2652   PLogEventEnd(MAT_IncompleteCholeskyFactorSymbolic,mat,perm,0,0);
2653   PetscFunctionReturn(0);
2654 }
2655 
2656 #undef __FUNC__
2657 #define __FUNC__ "MatGetArray"
2658 /*@C
2659    MatGetArray - Returns a pointer to the element values in the matrix.
2660    The result of this routine is dependent on the underlying matrix data-structure,
2661    and may not even work for certain matrix types. You MUST call MatRestoreArray() when you no
2662    longer need to access the array.
2663 
2664    Input Parameter:
2665 .  mat - the matrix
2666 
2667    Output Parameter:
2668 .  v - the location of the values
2669 
2670    Not Collective
2671 
2672    Currently only returns an array for the dense formats, giving access to the local portion
2673    of the matrix in the usual Fortran column oriented format.
2674 
2675   Fortran Note:
2676    The Fortran interface is slightly different from that given below.
2677    See the Fortran chapter of the users manual and
2678    petsc/src/mat/examples for details.
2679 
2680 .keywords: matrix, array, elements, values
2681 
2682 .seealso: MatRestoreArray()
2683 @*/
2684 int MatGetArray(Mat mat,Scalar **v)
2685 {
2686   int ierr;
2687 
2688   PetscFunctionBegin;
2689   PetscValidHeaderSpecific(mat,MAT_COOKIE);
2690   PetscValidPointer(v);
2691   if (!mat->ops->getarray) SETERRQ(PETSC_ERR_SUP,0,"");
2692   ierr = (*mat->ops->getarray)(mat,v);CHKERRQ(ierr);
2693   PetscFunctionReturn(0);
2694 }
2695 
2696 #undef __FUNC__
2697 #define __FUNC__ "MatRestoreArray"
2698 /*@C
2699    MatRestoreArray - Restores the matrix after MatGetArray() has been called.
2700 
2701    Input Parameter:
2702 .  mat - the matrix
2703 .  v - the location of the values
2704 
2705    Not Collective
2706 
2707     Fortran Note:
2708    The Fortran interface is slightly different from that given below.
2709    See the users manual and petsc/src/mat/examples for details.
2710 
2711 .keywords: matrix, array, elements, values, restore
2712 
2713 .seealso: MatGetArray()
2714 @*/
2715 int MatRestoreArray(Mat mat,Scalar **v)
2716 {
2717   int ierr;
2718 
2719   PetscFunctionBegin;
2720   PetscValidHeaderSpecific(mat,MAT_COOKIE);
2721   PetscValidPointer(v);
2722   if (!mat->ops->restorearray) SETERRQ(PETSC_ERR_SUP,0,"");
2723   ierr = (*mat->ops->restorearray)(mat,v);CHKERRQ(ierr);
2724   PetscFunctionReturn(0);
2725 }
2726 
2727 #undef __FUNC__
2728 #define __FUNC__ "MatGetSubMatrices"
2729 /*@C
2730    MatGetSubMatrices - Extracts several submatrices from a matrix. If submat
2731      points to an array of valid matrices, they may be reused to store the new
2732      submatrices.
2733 
2734    Input Parameters:
2735 .  mat - the matrix
2736 .  n   - the number of submatrixes to be extracted
2737 .  irow, icol - index sets of rows and columns to extract
2738 .  scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
2739 
2740    Output Parameter:
2741 .  submat - the array of submatrices
2742 
2743    Collective on Mat
2744 
2745    Notes:
2746    MatGetSubMatrices() can extract only sequential submatrices
2747    (from both sequential and parallel matrices). Use MatGetSubMatrix()
2748    to extract a parallel submatrix.
2749 
2750    When extracting submatrices from a parallel matrix, each processor can
2751    form a different submatrix by setting the rows and columns of its
2752    individual index sets according to the local submatrix desired.
2753 
2754    When finished using the submatrices, the user should destroy
2755    them with MatDestroySubMatrices().
2756 
2757 .keywords: matrix, get, submatrix, submatrices
2758 
2759 .seealso: MatDestroyMatrices(), MatGetSubMatrix()
2760 @*/
2761 int MatGetSubMatrices(Mat mat,int n,IS *irow,IS *icol,MatGetSubMatrixCall scall,Mat **submat)
2762 {
2763   int ierr;
2764 
2765   PetscFunctionBegin;
2766   PetscValidHeaderSpecific(mat,MAT_COOKIE);
2767   if (!mat->ops->getsubmatrices) SETERRQ(PETSC_ERR_SUP,0,"");
2768   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix");
2769 
2770   PLogEventBegin(MAT_GetSubMatrices,mat,0,0,0);
2771   ierr = (*mat->ops->getsubmatrices)(mat,n,irow,icol,scall,submat); CHKERRQ(ierr);
2772   PLogEventEnd(MAT_GetSubMatrices,mat,0,0,0);
2773 
2774   PetscFunctionReturn(0);
2775 }
2776 
2777 #undef __FUNC__
2778 #define __FUNC__ "MatDestroyMatrices"
2779 /*@C
2780    MatDestroyMatrices - Destroys a set of matrices obtained with MatGetSubMatrices().
2781 
2782    Input Parameters:
2783 .  n - the number of local matrices
2784 .  mat - the matrices
2785 
2786    Collective on Mat
2787 
2788 .keywords: matrix, destroy, submatrix, submatrices
2789 
2790 .seealso: MatGetSubMatrices()
2791 @*/
2792 int MatDestroyMatrices(int n,Mat **mat)
2793 {
2794   int ierr,i;
2795 
2796   PetscFunctionBegin;
2797   if (n < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,1,"Trying to destroy negative number of matrices");
2798   PetscValidPointer(mat);
2799   for ( i=0; i<n; i++ ) {
2800     ierr = MatDestroy((*mat)[i]); CHKERRQ(ierr);
2801   }
2802   if (n) PetscFree(*mat);
2803   PetscFunctionReturn(0);
2804 }
2805 
2806 #undef __FUNC__
2807 #define __FUNC__ "MatIncreaseOverlap"
2808 /*@
2809    MatIncreaseOverlap - Given a set of submatrices indicated by index sets,
2810    replaces the index sets by larger ones that represent submatrices with
2811    additional overlap.
2812 
2813    Input Parameters:
2814 .  mat - the matrix
2815 .  n   - the number of index sets
2816 .  is  - the array of pointers to index sets
2817 .  ov  - the additional overlap requested
2818 
2819    Collective on Mat
2820 
2821 .keywords: matrix, overlap, Schwarz
2822 
2823 .seealso: MatGetSubMatrices()
2824 @*/
2825 int MatIncreaseOverlap(Mat mat,int n, IS *is,int ov)
2826 {
2827   int ierr;
2828 
2829   PetscFunctionBegin;
2830   PetscValidHeaderSpecific(mat,MAT_COOKIE);
2831   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix");
2832   if (mat->factor)     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix");
2833 
2834   if (ov == 0) PetscFunctionReturn(0);
2835   if (!mat->ops->increaseoverlap) SETERRQ(PETSC_ERR_SUP,0,"");
2836   PLogEventBegin(MAT_IncreaseOverlap,mat,0,0,0);
2837   ierr = (*mat->ops->increaseoverlap)(mat,n,is,ov); CHKERRQ(ierr);
2838   PLogEventEnd(MAT_IncreaseOverlap,mat,0,0,0);
2839   PetscFunctionReturn(0);
2840 }
2841 
2842 #undef __FUNC__
2843 #define __FUNC__ "MatPrintHelp"
2844 /*@
2845    MatPrintHelp - Prints all the options for the matrix.
2846 
2847    Input Parameter:
2848 .  mat - the matrix
2849 
2850    Collective on Mat
2851 
2852    Options Database Keys:
2853 $  -help, -h
2854 
2855 .keywords: mat, help
2856 
2857 .seealso: MatCreate(), MatCreateXXX()
2858 @*/
2859 int MatPrintHelp(Mat mat)
2860 {
2861   static int called = 0;
2862   int        ierr;
2863   MPI_Comm   comm;
2864 
2865   PetscFunctionBegin;
2866   PetscValidHeaderSpecific(mat,MAT_COOKIE);
2867 
2868   comm = mat->comm;
2869   if (!called) {
2870     (*PetscHelpPrintf)(comm,"General matrix options:\n");
2871     (*PetscHelpPrintf)(comm,"  -mat_view_info: view basic matrix info during MatAssemblyEnd()\n");
2872     (*PetscHelpPrintf)(comm,"  -mat_view_info_detailed: view detailed matrix info during MatAssemblyEnd()\n");
2873     (*PetscHelpPrintf)(comm,"  -mat_view_draw: draw nonzero matrix structure during MatAssemblyEnd()\n");
2874     (*PetscHelpPrintf)(comm,"      -draw_pause <sec>: set seconds of display pause\n");
2875     (*PetscHelpPrintf)(comm,"      -display <name>: set alternate display\n");
2876     called = 1;
2877   }
2878   if (mat->ops->printhelp) {
2879     ierr = (*mat->ops->printhelp)(mat); CHKERRQ(ierr);
2880   }
2881   PetscFunctionReturn(0);
2882 }
2883 
2884 #undef __FUNC__
2885 #define __FUNC__ "MatGetBlockSize"
2886 /*@
2887    MatGetBlockSize - Returns the matrix block size; useful especially for the
2888    block row and block diagonal formats.
2889 
2890    Input Parameter:
2891 .  mat - the matrix
2892 
2893    Output Parameter:
2894 .  bs - block size
2895 
2896    Not Collective
2897 
2898    Notes:
2899 $  block diagonal formats: MATSEQBDIAG, MATMPIBDIAG
2900 $  block row formats: MATSEQBAIJ, MATMPIBAIJ
2901 
2902 .keywords: matrix, get, block, size
2903 
2904 .seealso: MatCreateSeqBAIJ(), MatCreateMPIBAIJ(), MatCreateSeqBDiag(), MatCreateMPIBDiag()
2905 @*/
2906 int MatGetBlockSize(Mat mat,int *bs)
2907 {
2908   int ierr;
2909 
2910   PetscFunctionBegin;
2911   PetscValidHeaderSpecific(mat,MAT_COOKIE);
2912   PetscValidIntPointer(bs);
2913   if (!mat->ops->getblocksize) SETERRQ(PETSC_ERR_SUP,0,"");
2914   ierr = (*mat->ops->getblocksize)(mat,bs);CHKERRQ(ierr);
2915   PetscFunctionReturn(0);
2916 }
2917 
2918 #undef __FUNC__
2919 #define __FUNC__ "MatGetRowIJ"
2920 /*@C
2921       MatGetRowIJ - Returns the compress row storage i and j indices for sequential matrices.
2922                  EXPERTS ONLY.
2923 
2924   Input Parameters:
2925 .   mat - the matrix
2926 .   shift -  0 or 1 indicating we want the indices starting at 0 or 1
2927 .   symmetric - PETSC_TRUE or PETSC_FALSE indicating the matrix data structure should be
2928                 symmetrized
2929 
2930   Output Parameters:
2931 .   n - number of rows and columns in the (possibly compressed) matrix
2932 .   ia - the row indices
2933 .   ja - the column indices
2934 .   done - PETSC_TRUE or PETSC_FALSE indicated that the values have been returned
2935 
2936    Collective on Mat
2937 
2938 @*/
2939 int MatGetRowIJ(Mat mat,int shift,PetscTruth symmetric,int *n,int **ia,int** ja,PetscTruth *done)
2940 {
2941   int ierr;
2942 
2943   PetscFunctionBegin;
2944   PetscValidHeaderSpecific(mat,MAT_COOKIE);
2945   if (ia) PetscValidIntPointer(ia);
2946   if (ja) PetscValidIntPointer(ja);
2947   PetscValidIntPointer(done);
2948   if (!mat->ops->getrowij) *done = PETSC_FALSE;
2949   else {
2950     *done = PETSC_TRUE;
2951     ierr  = (*mat->ops->getrowij)(mat,shift,symmetric,n,ia,ja,done); CHKERRQ(ierr);
2952   }
2953   PetscFunctionReturn(0);
2954 }
2955 
2956 #undef __FUNC__
2957 #define __FUNC__ "MatGetColumnIJ"
2958 /*@C
2959       MatGetColumnIJ - Returns the compress Column storage i and j indices for sequential matrices.
2960                  EXPERTS ONLY.
2961 
2962   Input Parameters:
2963 .   mat - the matrix
2964 .   shift - 1 or zero indicating we want the indices starting at 0 or 1
2965 .   symmetric - PETSC_TRUE or PETSC_FALSE indicating the matrix data structure should be
2966                 symmetrized
2967 
2968   Output Parameters:
2969 .   n - number of Columns and columns in the (possibly compressed) matrix
2970 .   ia - the Column indices
2971 .   ja - the column indices
2972 .   done - PETSC_TRUE or PETSC_FALSE indicated that the values have been returned
2973 
2974    Collective on Mat
2975 @*/
2976 int MatGetColumnIJ(Mat mat,int shift,PetscTruth symmetric,int *n,int **ia,int** ja,PetscTruth *done)
2977 {
2978   int ierr;
2979 
2980   PetscFunctionBegin;
2981   PetscValidHeaderSpecific(mat,MAT_COOKIE);
2982   if (ia) PetscValidIntPointer(ia);
2983   if (ja) PetscValidIntPointer(ja);
2984   PetscValidIntPointer(done);
2985 
2986   if (!mat->ops->getcolumnij) *done = PETSC_FALSE;
2987   else {
2988     *done = PETSC_TRUE;
2989     ierr  = (*mat->ops->getcolumnij)(mat,shift,symmetric,n,ia,ja,done); CHKERRQ(ierr);
2990   }
2991   PetscFunctionReturn(0);
2992 }
2993 
2994 #undef __FUNC__
2995 #define __FUNC__ "MatRestoreRowIJ"
2996 /*@C
2997       MatRestoreRowIJ - Call after you are completed with the ia,ja indices obtained with
2998                      MatGetRowIJ(). EXPERTS ONLY.
2999 
3000   Input Parameters:
3001 .   mat - the matrix
3002 .   shift - 1 or zero indicating we want the indices starting at 0 or 1
3003 .   symmetric - PETSC_TRUE or PETSC_FALSE indicating the matrix data structure should be
3004                 symmetrized
3005 
3006   Output Parameters:
3007 .   n - size of (possibly compressed) matrix
3008 .   ia - the row indices
3009 .   ja - the column indices
3010 .   done - PETSC_TRUE or PETSC_FALSE indicated that the values have been returned
3011 
3012   Collective on Mat
3013 
3014 @*/
3015 int MatRestoreRowIJ(Mat mat,int shift,PetscTruth symmetric,int *n,int **ia,int** ja,PetscTruth *done)
3016 {
3017   int ierr;
3018 
3019   PetscFunctionBegin;
3020   PetscValidHeaderSpecific(mat,MAT_COOKIE);
3021   if (ia) PetscValidIntPointer(ia);
3022   if (ja) PetscValidIntPointer(ja);
3023   PetscValidIntPointer(done);
3024 
3025   if (!mat->ops->restorerowij) *done = PETSC_FALSE;
3026   else {
3027     *done = PETSC_TRUE;
3028     ierr  = (*mat->ops->restorerowij)(mat,shift,symmetric,n,ia,ja,done); CHKERRQ(ierr);
3029   }
3030   PetscFunctionReturn(0);
3031 }
3032 
3033 #undef __FUNC__
3034 #define __FUNC__ "MatRestoreColumnIJ"
3035 /*@C
3036       MatRestoreColumnIJ - Call after you are completed with the ia,ja indices obtained with
3037                      MatGetColumnIJ(). EXPERTS ONLY.
3038 
3039   Input Parameters:
3040 .   mat - the matrix
3041 .   shift - 1 or zero indicating we want the indices starting at 0 or 1
3042 .   symmetric - PETSC_TRUE or PETSC_FALSE indicating the matrix data structure should be
3043                 symmetrized
3044 
3045   Output Parameters:
3046 .   n - size of (possibly compressed) matrix
3047 .   ia - the Column indices
3048 .   ja - the column indices
3049 .   done - PETSC_TRUE or PETSC_FALSE indicated that the values have been returned
3050 
3051   Collective on Mat
3052 
3053 @*/
3054 int MatRestoreColumnIJ(Mat mat,int shift,PetscTruth symmetric,int *n,int **ia,int** ja,PetscTruth *done)
3055 {
3056   int ierr;
3057 
3058   PetscFunctionBegin;
3059   PetscValidHeaderSpecific(mat,MAT_COOKIE);
3060   if (ia) PetscValidIntPointer(ia);
3061   if (ja) PetscValidIntPointer(ja);
3062   PetscValidIntPointer(done);
3063 
3064   if (!mat->ops->restorecolumnij) *done = PETSC_FALSE;
3065   else {
3066     *done = PETSC_TRUE;
3067     ierr  = (*mat->ops->restorecolumnij)(mat,shift,symmetric,n,ia,ja,done); CHKERRQ(ierr);
3068   }
3069   PetscFunctionReturn(0);
3070 }
3071 
3072 #undef __FUNC__
3073 #define __FUNC__ "MatColoringPatch"
3074 /*@C
3075       MatColoringPatch - EXPERTS ONLY, used inside matrix coloring routines that
3076           use matGetRowIJ() and/or MatGetColumnIJ().
3077 
3078   Input Parameters:
3079 .   mat - the matrix
3080 .   n   - number of colors
3081 .   colorarray - array indicating color for each column
3082 
3083   Output Parameters:
3084 .   iscoloring - coloring generated using colorarray information
3085 
3086   Collective on Mat
3087 
3088 @*/
3089 int MatColoringPatch(Mat mat,int n,int *colorarray,ISColoring *iscoloring)
3090 {
3091   int ierr;
3092 
3093   PetscFunctionBegin;
3094   PetscValidHeaderSpecific(mat,MAT_COOKIE);
3095   PetscValidIntPointer(colorarray);
3096 
3097   if (!mat->ops->coloringpatch) {SETERRQ(PETSC_ERR_SUP,0,"");}
3098   else {
3099     ierr  = (*mat->ops->coloringpatch)(mat,n,colorarray,iscoloring); CHKERRQ(ierr);
3100   }
3101   PetscFunctionReturn(0);
3102 }
3103 
3104 
3105 #undef __FUNC__
3106 #define __FUNC__ "MatSetUnfactored"
3107 /*@
3108    MatSetUnfactored - Resets a factored matrix to be treated as unfactored.
3109 
3110    Input Paramter:
3111 .  mat - the factored matrix to be reset
3112 
3113    Collective on Mat
3114 
3115    Notes:
3116    This routine should be used only with factored matrices formed by in-place
3117    factorization via ILU(0) (or by in-place LU factorization for the MATSEQDENSE
3118    format).  This option can save memory, for example, when solving nonlinear
3119    systems with a matrix-free Newton-Krylov method and a matrix-based, in-place
3120    ILU(0) preconditioner.
3121 
3122   Note that one can specify in-place ILU(0) factorization by calling
3123 $     PCType(pc,PCILU);
3124 $     PCILUSeUseInPlace(pc);
3125   or by using the options -pc_type ilu -pc_ilu_in_place
3126 
3127   In-place factorization ILU(0) can also be used as a local
3128   solver for the blocks within the block Jacobi or additive Schwarz
3129   methods (runtime option: -sub_pc_ilu_in_place).  See the discussion
3130   of these preconditioners in the users manual for details on setting
3131   local solver options.
3132 
3133 .seealso: PCILUSetUseInPlace(), PCLUSetUseInPlace()
3134 
3135 .keywords: matrix-free, in-place ILU, in-place LU
3136 @*/
3137 int MatSetUnfactored(Mat mat)
3138 {
3139   int ierr;
3140 
3141   PetscFunctionBegin;
3142   PetscValidHeaderSpecific(mat,MAT_COOKIE);
3143   mat->factor = 0;
3144   if (!mat->ops->setunfactored) PetscFunctionReturn(0);
3145   ierr = (*mat->ops->setunfactored)(mat); CHKERRQ(ierr);
3146   PetscFunctionReturn(0);
3147 }
3148 
3149 #undef __FUNC__
3150 #define __FUNC__ "MatGetType"
3151 /*@C
3152    MatGetType - Gets the matrix type and name (as a string) from the matrix.
3153 
3154    Input Parameter:
3155 .  mat - the matrix
3156 
3157    Output Parameter:
3158 .  type - the matrix type (or use PETSC_NULL)
3159 .  name - name of matrix type (or use PETSC_NULL)
3160 
3161    Not Collective
3162 
3163 .keywords: matrix, get, type, name
3164 @*/
3165 int MatGetType(Mat mat,MatType *type,char **name)
3166 {
3167   int  itype = (int)mat->type;
3168   char *matname[10];
3169 
3170   PetscFunctionBegin;
3171   PetscValidHeaderSpecific(mat,MAT_COOKIE);
3172 
3173   if (type) *type = (MatType) mat->type;
3174   if (name) {
3175     /* Note:  Be sure that this list corresponds to the enum in mat.h */
3176     matname[0] = "MATSEQDENSE";
3177     matname[1] = "MATSEQAIJ";
3178     matname[2] = "MATMPIAIJ";
3179     matname[3] = "MATSHELL";
3180     matname[4] = "MATMPIROWBS";
3181     matname[5] = "MATSEQBDIAG";
3182     matname[6] = "MATMPIBDIAG";
3183     matname[7] = "MATMPIDENSE";
3184     matname[8] = "MATSEQBAIJ";
3185     matname[9] = "MATMPIBAIJ";
3186 
3187     if (itype < 0 || itype > 9) *name = "Unknown matrix type";
3188     else                        *name = matname[itype];
3189   }
3190   PetscFunctionReturn(0);
3191 }
3192 
3193 /*MC
3194     MatGetArrayF90 - Accesses a matrix array from Fortran90.
3195 
3196     Input Parameter:
3197 .   x - matrix
3198 
3199     Output Parameter:
3200 .   xx_v - the Fortran90 pointer to the array
3201 .   ierr - error code
3202 
3203     Synopsis:
3204     MatGetArrayF90(Mat x,{Scalar, pointer :: xx_v(:)},integer ierr)
3205 
3206     Usage:
3207 $     Scalar, pointer xx_v(:)
3208 $     ....
3209 $     call MatGetArrayF90(x,xx_v,ierr)
3210 $     a = xx_v(3)
3211 $     call MatRestoreArrayF90(x,xx_v,ierr)
3212 
3213     Notes:
3214     Currently only supported using the NAG F90 compiler.
3215 
3216 .seealso:  MatRestoreArrayF90(), MatGetArray(), MatRestoreArray()
3217 
3218 .keywords:  matrix, array, f90
3219 M*/
3220 
3221 /*MC
3222     MatRestoreArrayF90 - Restores a matrix array that has been
3223     accessed with MatGetArrayF90().
3224 
3225     Input Parameters:
3226 .   x - matrix
3227 .   xx_v - the Fortran90 pointer to the array
3228 
3229     Output Parameter:
3230 .   ierr - error code
3231 
3232     Synopsis:
3233     MatRestoreArrayF90(Mat x,{Scalar, pointer :: xx_v(:)},integer ierr)
3234 
3235     Example of Usage:
3236 $      Scalar, pointer xx_v(:)
3237 $      ....
3238 $      call MatGetArrayF90(x,xx_v,ierr)
3239 $      a = xx_v(3)
3240 $      call MatRestoreArrayF90(x,xx_v,ierr)
3241 
3242     Notes:
3243     Currently only supported using the NAG F90 compiler.
3244 
3245 .seealso:  MatGetArrayF90(), MatGetArray(), MatRestoreArray()
3246 
3247 .keywords:  matrix, array, f90
3248 M*/
3249 
3250 
3251 #undef __FUNC__
3252 #define __FUNC__ "MatGetSubMatrix"
3253 /*@
3254      MatGetSubMatrix - Gets a single submatrix on the same number of processors
3255                        as the original matrix.
3256 
3257    Input Parameters:
3258 .   mat - the original matrix
3259 .   isrow - rows this processor should obtain
3260 .   iscol - columns for all processors you wish kept
3261 .   csize - number of columns "local" to this processor (does nothing for sequential
3262             matrices). This should match the result from VecGetLocalSize() if you
3263             plan to use the matrix in a A*x
3264 .  cll - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
3265 
3266    Output Parameters:
3267 .   newmat - the new submatrix, of the same type as the old
3268 
3269    Collective on Mat
3270 
3271 .seealso: MatGetSubMatrices()
3272 
3273 @*/
3274 int MatGetSubMatrix(Mat mat,IS isrow,IS iscol,int csize,MatGetSubMatrixCall cll,Mat *newmat)
3275 {
3276   int     ierr, size;
3277   Mat     *local;
3278 
3279   PetscFunctionBegin;
3280   MPI_Comm_size(mat->comm,&size);
3281 
3282   /* if original matrix is on just one processor then use submatrix generated */
3283   if (!mat->ops->getsubmatrix && size == 1 && cll == MAT_REUSE_MATRIX) {
3284     ierr = MatGetSubMatrices(mat,1,&isrow,&iscol,MAT_REUSE_MATRIX,&newmat);CHKERRQ(ierr);
3285     PetscFunctionReturn(0);
3286   } else if (!mat->ops->getsubmatrix && size == 1) {
3287     ierr = MatGetSubMatrices(mat,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&local);CHKERRQ(ierr);
3288     *newmat = *local;
3289     PetscFree(local);
3290     PetscFunctionReturn(0);
3291   }
3292 
3293   if (!mat->ops->getsubmatrix) SETERRQ(PETSC_ERR_SUP,0,"Not currently implemented");
3294   ierr = (*mat->ops->getsubmatrix)(mat,isrow,iscol,csize,cll,newmat);CHKERRQ(ierr);
3295   PetscFunctionReturn(0);
3296 }
3297 
3298 
3299 
3300