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