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