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