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