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