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