xref: /petsc/src/mat/interface/matrix.c (revision 60625e5193e1faa2b4b90582eb81189a9b39b613)
1 /*$Id: matrix.c,v 1.394 2001/03/09 18:59:52 balay Exp bsmith $*/
2 
3 /*
4    This is where the abstract matrix operations are defined
5 */
6 
7 #include "src/mat/matimpl.h"        /*I "petscmat.h" I*/
8 #include "src/vec/vecimpl.h"
9 
10 #undef __FUNC__
11 #define __FUNC__ "MatGetRow"
12 /*@C
13    MatGetRow - Gets a row of a matrix.  You MUST call MatRestoreRow()
14    for each row that you get to ensure that your application does
15    not bleed memory.
16 
17    Not Collective
18 
19    Input Parameters:
20 +  mat - the matrix
21 -  row - the row to get
22 
23    Output Parameters:
24 +  ncols -  the number of nonzeros in the row
25 .  cols - if not NULL, the column numbers
26 -  vals - if not NULL, the values
27 
28    Notes:
29    This routine is provided for people who need to have direct access
30    to the structure of a matrix.  We hope that we provide enough
31    high-level matrix routines that few users will need it.
32 
33    MatGetRow() always returns 0-based column indices, regardless of
34    whether the internal representation is 0-based (default) or 1-based.
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    You can only have one call to MatGetRow() outstanding for a particular
44    matrix at a time, per processor. MatGetRow() can only obtained rows
45    associated with the given processor, it cannot get rows from the
46    other processors; for that we suggest using MatGetSubMatrices(), then
47    MatGetRow() on the submatrix. The row indix passed to MatGetRows()
48    is in the global number of rows.
49 
50    Fortran Notes:
51    The calling sequence from Fortran is
52 .vb
53    MatGetRow(matrix,row,ncols,cols,values,ierr)
54          Mat     matrix (input)
55          integer row    (input)
56          integer ncols  (output)
57          integer cols(maxcols) (output)
58          double precision (or double complex) values(maxcols) output
59 .ve
60    where maxcols >= maximum nonzeros in any row of the matrix.
61 
62    Caution:
63    Do not try to change the contents of the output arrays (cols and vals).
64    In some cases, this may corrupt the matrix.
65 
66    Level: advanced
67 
68    Concepts: matrices^row access
69 
70 .seealso: MatRestoreRow(), MatSetValues(), MatGetValues(), MatGetSubmatrices(), MatGetDiagonal()
71 @*/
72 int MatGetRow(Mat mat,int row,int *ncols,int **cols,Scalar **vals)
73 {
74   int   ierr;
75 
76   PetscFunctionBegin;
77   PetscValidHeaderSpecific(mat,MAT_COOKIE);
78   PetscValidIntPointer(ncols);
79   PetscValidType(mat);
80   MatPreallocated(mat);
81   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
82   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
83   if (!mat->ops->getrow) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name);
84   ierr = PetscLogEventBegin(MAT_GetRow,mat,0,0,0);CHKERRQ(ierr);
85   ierr = (*mat->ops->getrow)(mat,row,ncols,cols,vals);CHKERRQ(ierr);
86   ierr = PetscLogEventEnd(MAT_GetRow,mat,0,0,0);CHKERRQ(ierr);
87   PetscFunctionReturn(0);
88 }
89 
90 #undef __FUNC__
91 #define __FUNC__ "MatRestoreRow"
92 /*@C
93    MatRestoreRow - Frees any temporary space allocated by MatGetRow().
94 
95    Not Collective
96 
97    Input Parameters:
98 +  mat - the matrix
99 .  row - the row to get
100 .  ncols, cols - the number of nonzeros and their columns
101 -  vals - if nonzero the column values
102 
103    Notes:
104    This routine should be called after you have finished examining the entries.
105 
106    Fortran Notes:
107    The calling sequence from Fortran is
108 .vb
109    MatRestoreRow(matrix,row,ncols,cols,values,ierr)
110       Mat     matrix (input)
111       integer row    (input)
112       integer ncols  (output)
113       integer cols(maxcols) (output)
114       double precision (or double complex) values(maxcols) output
115 .ve
116    Where maxcols >= maximum nonzeros in any row of the matrix.
117 
118    In Fortran MatRestoreRow() MUST be called after MatGetRow()
119    before another call to MatGetRow() can be made.
120 
121    Level: advanced
122 
123 .seealso:  MatGetRow()
124 @*/
125 int MatRestoreRow(Mat mat,int row,int *ncols,int **cols,Scalar **vals)
126 {
127   int ierr;
128 
129   PetscFunctionBegin;
130   PetscValidHeaderSpecific(mat,MAT_COOKIE);
131   PetscValidIntPointer(ncols);
132   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
133   if (!mat->ops->restorerow) PetscFunctionReturn(0);
134   ierr = (*mat->ops->restorerow)(mat,row,ncols,cols,vals);CHKERRQ(ierr);
135   PetscFunctionReturn(0);
136 }
137 
138 #undef __FUNC__
139 #define __FUNC__ "MatView"
140 /*@C
141    MatView - Visualizes a matrix object.
142 
143    Collective on Mat
144 
145    Input Parameters:
146 +  mat - the matrix
147 -  ptr - visualization context
148 
149   Notes:
150   The available visualization contexts include
151 +    PETSC_VIEWER_STDOUT_SELF - standard output (default)
152 .    PETSC_VIEWER_STDOUT_WORLD - synchronized standard
153         output where only the first processor opens
154         the file.  All other processors send their
155         data to the first processor to print.
156 -     PETSC_VIEWER_DRAW_WORLD - graphical display of nonzero structure
157 
158    The user can open alternative visualization contexts with
159 +    PetscViewerASCIIOpen() - Outputs matrix to a specified file
160 .    PetscViewerBinaryOpen() - Outputs matrix in binary to a
161          specified file; corresponding input uses MatLoad()
162 .    PetscViewerDrawOpen() - Outputs nonzero matrix structure to
163          an X window display
164 -    PetscViewerSocketOpen() - Outputs matrix to Socket viewer.
165          Currently only the sequential dense and AIJ
166          matrix types support the Socket viewer.
167 
168    The user can call PetscViewerSetFormat() to specify the output
169    format of ASCII printed objects (when using PETSC_VIEWER_STDOUT_SELF,
170    PETSC_VIEWER_STDOUT_WORLD and PetscViewerASCIIOpen).  Available formats include
171 +    PETSC_VIEWER_ASCII_DEFAULT - default, prints matrix contents
172 .    PETSC_VIEWER_ASCII_MATLAB - prints matrix contents in Matlab format
173 .    PETSC_VIEWER_ASCII_DENSE - prints entire matrix including zeros
174 .    PETSC_VIEWER_ASCII_COMMON - prints matrix contents, using a sparse
175          format common among all matrix types
176 .    PETSC_VIEWER_ASCII_IMPL - prints matrix contents, using an implementation-specific
177          format (which is in many cases the same as the default)
178 .    PETSC_VIEWER_ASCII_INFO - prints basic information about the matrix
179          size and structure (not the matrix entries)
180 -    PETSC_VIEWER_ASCII_INFO_LONG - prints more detailed information about
181          the matrix structure
182 
183    Level: beginner
184 
185    Concepts: matrices^viewing
186    Concepts: matrices^plotting
187    Concepts: matrices^printing
188 
189 .seealso: PetscViewerSetFormat(), PetscViewerASCIIOpen(), PetscViewerDrawOpen(),
190           PetscViewerSocketOpen(), PetscViewerBinaryOpen(), MatLoad()
191 @*/
192 int MatView(Mat mat,PetscViewer viewer)
193 {
194   int               ierr,rows,cols;
195   PetscTruth        isascii;
196   char              *cstr;
197   PetscViewerFormat format;
198 
199   PetscFunctionBegin;
200   PetscValidHeaderSpecific(mat,MAT_COOKIE);
201   PetscValidType(mat);
202   MatPreallocated(mat);
203   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(mat->comm);
204   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE);
205   PetscCheckSameComm(mat,viewer);
206   if (!mat->assembled) SETERRQ(1,"Must call MatAssemblyBegin/End() before viewing matrix");
207 
208   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
209   if (isascii) {
210     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
211     if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_LONG) {
212       ierr = PetscViewerASCIIPrintf(viewer,"Matrix Object:\n");CHKERRQ(ierr);
213       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
214       ierr = MatGetType(mat,&cstr);CHKERRQ(ierr);
215       ierr = MatGetSize(mat,&rows,&cols);CHKERRQ(ierr);
216       ierr = PetscViewerASCIIPrintf(viewer,"type=%s, rows=%d, cols=%d\n",cstr,rows,cols);CHKERRQ(ierr);
217       if (mat->ops->getinfo) {
218         MatInfo info;
219         ierr = MatGetInfo(mat,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr);
220         ierr = PetscViewerASCIIPrintf(viewer,"total: nonzeros=%d, allocated nonzeros=%d\n",
221                           (int)info.nz_used,(int)info.nz_allocated);CHKERRQ(ierr);
222       }
223     }
224   }
225   if (mat->ops->view) {
226     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
227     ierr = (*mat->ops->view)(mat,viewer);CHKERRQ(ierr);
228     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
229   } else if (!isascii) {
230     SETERRQ1(1,"Viewer type %s not supported",((PetscObject)viewer)->type_name);
231   }
232   if (isascii) {
233     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
234     if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_LONG) {
235       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
236     }
237   }
238   PetscFunctionReturn(0);
239 }
240 
241 #undef __FUNC__
242 #define __FUNC__ "MatScaleSystem"
243 /*@C
244    MatScaleSystem - Scale a vector solution and right hand side to
245    match the scaling of a scaled matrix.
246 
247    Collective on Mat
248 
249    Input Parameter:
250 +  mat - the matrix
251 .  x - solution vector (or PETSC_NULL)
252 +  b - right hand side vector (or PETSC_NULL)
253 
254 
255    Notes:
256    For AIJ, BAIJ, and BDiag matrix formats, the matrices are not
257    internally scaled, so this does nothing. For MPIROWBS it
258    permutes and diagonally scales.
259 
260    The SLES methods automatically call this routine when required
261    (via PCPreSolve()) so it is rarely used directly.
262 
263    Level: Developer
264 
265    Concepts: matrices^scaling
266 
267 .seealso: MatUseScaledForm(), MatUnScaleSystem()
268 @*/
269 int MatScaleSystem(Mat mat,Vec x,Vec b)
270 {
271   int ierr;
272 
273   PetscFunctionBegin;
274   PetscValidHeaderSpecific(mat,MAT_COOKIE);
275   PetscValidType(mat);
276   MatPreallocated(mat);
277   if (x) {PetscValidHeaderSpecific(x,VEC_COOKIE);PetscCheckSameComm(mat,x);}
278   if (b) {PetscValidHeaderSpecific(b,VEC_COOKIE);PetscCheckSameComm(mat,b);}
279 
280   if (mat->ops->scalesystem) {
281     ierr = (*mat->ops->scalesystem)(mat,x,b);CHKERRQ(ierr);
282   }
283   PetscFunctionReturn(0);
284 }
285 
286 #undef __FUNC__
287 #define __FUNC__ "MatUnScaleSystem"
288 /*@C
289    MatUnScaleSystem - Unscales a vector solution and right hand side to
290    match the original scaling of a scaled matrix.
291 
292    Collective on Mat
293 
294    Input Parameter:
295 +  mat - the matrix
296 .  x - solution vector (or PETSC_NULL)
297 +  b - right hand side vector (or PETSC_NULL)
298 
299 
300    Notes:
301    For AIJ, BAIJ, and BDiag matrix formats, the matrices are not
302    internally scaled, so this does nothing. For MPIROWBS it
303    permutes and diagonally scales.
304 
305    The SLES methods automatically call this routine when required
306    (via PCPreSolve()) so it is rarely used directly.
307 
308    Level: Developer
309 
310 .seealso: MatUseScaledForm(), MatScaleSystem()
311 @*/
312 int MatUnScaleSystem(Mat mat,Vec x,Vec b)
313 {
314   int ierr;
315 
316   PetscFunctionBegin;
317   PetscValidHeaderSpecific(mat,MAT_COOKIE);
318   PetscValidType(mat);
319   MatPreallocated(mat);
320   if (x) {PetscValidHeaderSpecific(x,VEC_COOKIE);PetscCheckSameComm(mat,x);}
321   if (b) {PetscValidHeaderSpecific(b,VEC_COOKIE);PetscCheckSameComm(mat,b);}
322   if (mat->ops->unscalesystem) {
323     ierr = (*mat->ops->unscalesystem)(mat,x,b);CHKERRQ(ierr);
324   }
325   PetscFunctionReturn(0);
326 }
327 
328 #undef __FUNC__
329 #define __FUNC__ "MatUseScaledForm"
330 /*@C
331    MatUseScaledForm - For matrix storage formats that scale the
332    matrix (for example MPIRowBS matrices are diagonally scaled on
333    assembly) indicates matrix operations (MatMult() etc) are
334    applied using the scaled matrix.
335 
336    Collective on Mat
337 
338    Input Parameter:
339 +  mat - the matrix
340 -  scaled - PETSC_TRUE for applying the scaled, PETSC_FALSE for
341             applying the original matrix
342 
343    Notes:
344    For scaled matrix formats, applying the original, unscaled matrix
345    will be slightly more expensive
346 
347    Level: Developer
348 
349 .seealso: MatScaleSystem(), MatUnScaleSystem()
350 @*/
351 int MatUseScaledForm(Mat mat,PetscTruth scaled)
352 {
353   int ierr;
354 
355   PetscFunctionBegin;
356   PetscValidHeaderSpecific(mat,MAT_COOKIE);
357   PetscValidType(mat);
358   MatPreallocated(mat);
359   if (mat->ops->usescaledform) {
360     ierr = (*mat->ops->usescaledform)(mat,scaled);CHKERRQ(ierr);
361   }
362   PetscFunctionReturn(0);
363 }
364 
365 #undef __FUNC__
366 #define __FUNC__ "MatDestroy"
367 /*@C
368    MatDestroy - Frees space taken by a matrix.
369 
370    Collective on Mat
371 
372    Input Parameter:
373 .  A - the matrix
374 
375    Level: beginner
376 
377 @*/
378 int MatDestroy(Mat A)
379 {
380   int ierr;
381 
382   PetscFunctionBegin;
383   PetscValidHeaderSpecific(A,MAT_COOKIE);
384   PetscValidType(A);
385   MatPreallocated(A);
386   if (--A->refct > 0) PetscFunctionReturn(0);
387 
388   /* if memory was published with AMS then destroy it */
389   ierr = PetscObjectDepublish(A);CHKERRQ(ierr);
390   if (A->mapping) {
391     ierr = ISLocalToGlobalMappingDestroy(A->mapping);CHKERRQ(ierr);
392   }
393   if (A->bmapping) {
394     ierr = ISLocalToGlobalMappingDestroy(A->bmapping);CHKERRQ(ierr);
395   }
396   if (A->rmap) {
397     ierr = MapDestroy(A->rmap);CHKERRQ(ierr);
398   }
399   if (A->cmap) {
400     ierr = MapDestroy(A->cmap);CHKERRQ(ierr);
401   }
402 
403   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
404   PetscLogObjectDestroy(A);
405   PetscHeaderDestroy(A);
406   PetscFunctionReturn(0);
407 }
408 
409 #undef __FUNC__
410 #define __FUNC__ "MatValid"
411 /*@
412    MatValid - Checks whether a matrix object is valid.
413 
414    Collective on Mat
415 
416    Input Parameter:
417 .  m - the matrix to check
418 
419    Output Parameter:
420    flg - flag indicating matrix status, either
421    PETSC_TRUE if matrix is valid, or PETSC_FALSE otherwise.
422 
423    Level: developer
424 
425    Concepts: matrices^validity
426 @*/
427 int MatValid(Mat m,PetscTruth *flg)
428 {
429   PetscFunctionBegin;
430   PetscValidIntPointer(flg);
431   if (!m)                           *flg = PETSC_FALSE;
432   else if (m->cookie != MAT_COOKIE) *flg = PETSC_FALSE;
433   else                              *flg = PETSC_TRUE;
434   PetscFunctionReturn(0);
435 }
436 
437 #undef __FUNC__
438 #define __FUNC__ "MatSetValues"
439 /*@
440    MatSetValues - Inserts or adds a block of values into a matrix.
441    These values may be cached, so MatAssemblyBegin() and MatAssemblyEnd()
442    MUST be called after all calls to MatSetValues() have been completed.
443 
444    Not Collective
445 
446    Input Parameters:
447 +  mat - the matrix
448 .  v - a logically two-dimensional array of values
449 .  m, idxm - the number of rows and their global indices
450 .  n, idxn - the number of columns and their global indices
451 -  addv - either ADD_VALUES or INSERT_VALUES, where
452    ADD_VALUES adds values to any existing entries, and
453    INSERT_VALUES replaces existing entries with new values
454 
455    Notes:
456    By default the values, v, are row-oriented and unsorted.
457    See MatSetOption() for other options.
458 
459    Calls to MatSetValues() with the INSERT_VALUES and ADD_VALUES
460    options cannot be mixed without intervening calls to the assembly
461    routines.
462 
463    MatSetValues() uses 0-based row and column numbers in Fortran
464    as well as in C.
465 
466    Negative indices may be passed in idxm and idxn, these rows and columns are
467    simply ignored. This allows easily inserting element stiffness matrices
468    with homogeneous Dirchlet boundary conditions that you don't want represented
469    in the matrix.
470 
471    Efficiency Alert:
472    The routine MatSetValuesBlocked() may offer much better efficiency
473    for users of block sparse formats (MATSEQBAIJ and MATMPIBAIJ).
474 
475    Level: beginner
476 
477    Concepts: matrices^putting entries in
478 
479 .seealso: MatSetOption(), MatAssemblyBegin(), MatAssemblyEnd(), MatSetValuesBlocked(), MatSetValuesLocal()
480 @*/
481 int MatSetValues(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v,InsertMode addv)
482 {
483   int ierr;
484 
485   PetscFunctionBegin;
486   if (!m || !n) PetscFunctionReturn(0); /* no values to insert */
487   PetscValidHeaderSpecific(mat,MAT_COOKIE);
488   PetscValidType(mat);
489   MatPreallocated(mat);
490   PetscValidIntPointer(idxm);
491   PetscValidIntPointer(idxn);
492   PetscValidScalarPointer(v);
493   if (mat->insertmode == NOT_SET_VALUES) {
494     mat->insertmode = addv;
495   }
496 #if defined(PETSC_USE_BOPT_g)
497   else if (mat->insertmode != addv) {
498     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot mix add values and insert values");
499   }
500   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
501 #endif
502 
503   if (mat->assembled) {
504     mat->was_assembled = PETSC_TRUE;
505     mat->assembled     = PETSC_FALSE;
506   }
507   ierr = PetscLogEventBegin(MAT_SetValues,mat,0,0,0);CHKERRQ(ierr);
508   if (!mat->ops->setvalues) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name);
509   ierr = (*mat->ops->setvalues)(mat,m,idxm,n,idxn,v,addv);CHKERRQ(ierr);
510   ierr = PetscLogEventEnd(MAT_SetValues,mat,0,0,0);CHKERRQ(ierr);
511   PetscFunctionReturn(0);
512 }
513 
514 #undef __FUNC__
515 #define __FUNC__ "MatSetValuesStencil"
516 /*@
517    MatSetValuesStencil - Inserts or adds a block of values into a matrix.
518      Using structured grid indexing
519 
520    Not Collective
521 
522    Input Parameters:
523 +  mat - the matrix
524 .  v - a logically two-dimensional array of values
525 .  m - number of rows being entered
526 .  idxm - see below
527 .  n - number of columns being entered
528 .  idxn -  see below
529 -  addv - either ADD_VALUES or INSERT_VALUES, where
530    ADD_VALUES adds values to any existing entries, and
531    INSERT_VALUES replaces existing entries with new values
532 
533    Notes:
534    Format of idxm and idxn:
535 $    If dof > 1 let dim = dimension of the grid (1,2, or 3)+1;
536 $     then idxm[dim*i+dim-1] = component number
537 $        idxm[dim*i+dim-2]   = logical coordinate in x
538 $        idxm[dim*i+dim-3]   = logical coordinate in y (if it exists)
539 $        idxm[dim*i+dim-4]   = logical coordinate in z (if it exists)
540 $     else let dim = dimension of the grid (1,2, or 3);
541 $        idxm[dim*i+dim-1]   = logical coordinate in x
542 $        idxm[dim*i+dim-2]   = logical coordinate in y (if it exists)
543 $        idxm[dim*i+dim-3]   = logical coordinate in z (if it exists)
544 
545 
546    By default the values, v, are row-oriented and unsorted.
547    See MatSetOption() for other options.
548 
549    Calls to MatSetValuesStencil() with the INSERT_VALUES and ADD_VALUES
550    options cannot be mixed without intervening calls to the assembly
551    routines.
552 
553    MatSetValuesStencil() uses 0-based row and column numbers in Fortran
554    as well as in C.
555 
556    Negative indices may be passed in idxm and idxn, these rows and columns are
557    simply ignored. This allows easily inserting element stiffness matrices
558    with homogeneous Dirchlet boundary conditions that you don't want represented
559    in the matrix.
560 
561    Inspired by the structured grid interface to the HYPRE package
562    (www.llnl.gov/CASC/hyper)
563 
564    Efficiency Alert:
565    The routine MatSetValuesBlockedStencil() may offer much better efficiency
566    for users of block sparse formats (MATSEQBAIJ and MATMPIBAIJ).
567 
568    Level: beginner
569 
570    Concepts: matrices^putting entries in
571 
572 .seealso: MatSetOption(), MatAssemblyBegin(), MatAssemblyEnd(), MatSetValuesBlocked(), MatSetValuesLocal()
573           MatSetValues(), MatSetValuesBlockedStencil(), MatSetStencil()
574 @*/
575 int MatSetValuesStencil(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v,InsertMode addv)
576 {
577   int j,i,ierr,jdxm[128],jdxn[128],dim = mat->stencil.dim,*dims = mat->stencil.dims+1,tmp;
578 
579   PetscFunctionBegin;
580   if (!m || !n) PetscFunctionReturn(0); /* no values to insert */
581   PetscValidHeaderSpecific(mat,MAT_COOKIE);
582   PetscValidType(mat);
583   PetscValidIntPointer(idxm);
584   PetscValidIntPointer(idxn);
585   PetscValidScalarPointer(v);
586 
587   for (i=0; i<m; i++) {
588     tmp = idxm[dim*i];
589     for (j=0; j<dim-1; j++) {
590       tmp = tmp*dims[j] + idxm[dim*i+j+1];
591     }
592     jdxm[i] = tmp;
593   }
594   for (i=0; i<n; i++) {
595     tmp = idxn[dim*i];
596     for (j=0; j<dim-1; j++) {
597       tmp = tmp*dims[j] + idxn[dim*i+j+1];
598     }
599     jdxn[i] = tmp;
600   }
601   ierr = MatSetValues(mat,m,jdxm,n,jdxn,v,addv);CHKERRQ(ierr);
602   PetscFunctionReturn(0);
603 }
604 
605 #undef __FUNC__
606 #define __FUNC__ "MatSetValuesBlocked"
607 /*@
608    MatSetValuesBlocked - Inserts or adds a block of values into a matrix.
609 
610    Not Collective
611 
612    Input Parameters:
613 +  mat - the matrix
614 .  v - a logically two-dimensional array of values
615 .  m, idxm - the number of block rows and their global block indices
616 .  n, idxn - the number of block columns and their global block indices
617 -  addv - either ADD_VALUES or INSERT_VALUES, where
618    ADD_VALUES adds values to any existing entries, and
619    INSERT_VALUES replaces existing entries with new values
620 
621    Notes:
622    By default the values, v, are row-oriented and unsorted. So the layout of
623    v is the same as for MatSetValues(). See MatSetOption() for other options.
624 
625    Calls to MatSetValuesBlocked() with the INSERT_VALUES and ADD_VALUES
626    options cannot be mixed without intervening calls to the assembly
627    routines.
628 
629    MatSetValuesBlocked() uses 0-based row and column numbers in Fortran
630    as well as in C.
631 
632    Negative indices may be passed in idxm and idxn, these rows and columns are
633    simply ignored. This allows easily inserting element stiffness matrices
634    with homogeneous Dirchlet boundary conditions that you don't want represented
635    in the matrix.
636 
637    Each time an entry is set within a sparse matrix via MatSetValues(),
638    internal searching must be done to determine where to place the the
639    data in the matrix storage space.  By instead inserting blocks of
640    entries via MatSetValuesBlocked(), the overhead of matrix assembly is
641    reduced.
642 
643    Restrictions:
644    MatSetValuesBlocked() is currently supported only for the block AIJ
645    matrix format (MATSEQBAIJ and MATMPIBAIJ, which are created via
646    MatCreateSeqBAIJ() and MatCreateMPIBAIJ()).
647 
648    Level: intermediate
649 
650    Concepts: matrices^putting entries in blocked
651 
652 .seealso: MatSetOption(), MatAssemblyBegin(), MatAssemblyEnd(), MatSetValues(), MatSetValuesBlockedLocal()
653 @*/
654 int MatSetValuesBlocked(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v,InsertMode addv)
655 {
656   int ierr;
657 
658   PetscFunctionBegin;
659   if (!m || !n) PetscFunctionReturn(0); /* no values to insert */
660   PetscValidHeaderSpecific(mat,MAT_COOKIE);
661   PetscValidType(mat);
662   MatPreallocated(mat);
663   PetscValidIntPointer(idxm);
664   PetscValidIntPointer(idxn);
665   PetscValidScalarPointer(v);
666   if (mat->insertmode == NOT_SET_VALUES) {
667     mat->insertmode = addv;
668   }
669 #if defined(PETSC_USE_BOPT_g)
670   else if (mat->insertmode != addv) {
671     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot mix add values and insert values");
672   }
673   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
674 #endif
675 
676   if (mat->assembled) {
677     mat->was_assembled = PETSC_TRUE;
678     mat->assembled     = PETSC_FALSE;
679   }
680   ierr = PetscLogEventBegin(MAT_SetValues,mat,0,0,0);CHKERRQ(ierr);
681   if (!mat->ops->setvaluesblocked) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name);
682   ierr = (*mat->ops->setvaluesblocked)(mat,m,idxm,n,idxn,v,addv);CHKERRQ(ierr);
683   ierr = PetscLogEventEnd(MAT_SetValues,mat,0,0,0);CHKERRQ(ierr);
684   PetscFunctionReturn(0);
685 }
686 
687 /*MC
688    MatSetValue - Set a single entry into a matrix.
689 
690    Synopsis:
691    void MatSetValue(Mat m,int row,int col,Scalar value,InsertMode mode);
692 
693    Not collective
694 
695    Input Parameters:
696 +  m - the matrix
697 .  row - the row location of the entry
698 .  col - the column location of the entry
699 .  value - the value to insert
700 -  mode - either INSERT_VALUES or ADD_VALUES
701 
702    Notes:
703    For efficiency one should use MatSetValues() and set several or many
704    values simultaneously if possible.
705 
706    Note that VecSetValue() does NOT return an error code (since this
707    is checked internally).
708 
709    Level: beginner
710 
711 .seealso: MatSetValues()
712 M*/
713 
714 #undef __FUNC__
715 #define __FUNC__ "MatGetValues"
716 /*@
717    MatGetValues - Gets a block of values from a matrix.
718 
719    Not Collective; currently only returns a local block
720 
721    Input Parameters:
722 +  mat - the matrix
723 .  v - a logically two-dimensional array for storing the values
724 .  m, idxm - the number of rows and their global indices
725 -  n, idxn - the number of columns and their global indices
726 
727    Notes:
728    The user must allocate space (m*n Scalars) for the values, v.
729    The values, v, are then returned in a row-oriented format,
730    analogous to that used by default in MatSetValues().
731 
732    MatGetValues() uses 0-based row and column numbers in
733    Fortran as well as in C.
734 
735    MatGetValues() requires that the matrix has been assembled
736    with MatAssemblyBegin()/MatAssemblyEnd().  Thus, calls to
737    MatSetValues() and MatGetValues() CANNOT be made in succession
738    without intermediate matrix assembly.
739 
740    Level: advanced
741 
742    Concepts: matrices^accessing values
743 
744 .seealso: MatGetRow(), MatGetSubmatrices(), MatSetValues()
745 @*/
746 int MatGetValues(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v)
747 {
748   int ierr;
749 
750   PetscFunctionBegin;
751   PetscValidHeaderSpecific(mat,MAT_COOKIE);
752   PetscValidType(mat);
753   MatPreallocated(mat);
754   PetscValidIntPointer(idxm);
755   PetscValidIntPointer(idxn);
756   PetscValidScalarPointer(v);
757   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
758   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
759   if (!mat->ops->getvalues) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name);
760 
761   ierr = PetscLogEventBegin(MAT_GetValues,mat,0,0,0);CHKERRQ(ierr);
762   ierr = (*mat->ops->getvalues)(mat,m,idxm,n,idxn,v);CHKERRQ(ierr);
763   ierr = PetscLogEventEnd(MAT_GetValues,mat,0,0,0);CHKERRQ(ierr);
764   PetscFunctionReturn(0);
765 }
766 
767 #undef __FUNC__
768 #define __FUNC__ "MatSetLocalToGlobalMapping"
769 /*@
770    MatSetLocalToGlobalMapping - Sets a local-to-global numbering for use by
771    the routine MatSetValuesLocal() to allow users to insert matrix entries
772    using a local (per-processor) numbering.
773 
774    Not Collective
775 
776    Input Parameters:
777 +  x - the matrix
778 -  mapping - mapping created with ISLocalToGlobalMappingCreate()
779              or ISLocalToGlobalMappingCreateIS()
780 
781    Level: intermediate
782 
783    Concepts: matrices^local to global mapping
784    Concepts: local to global mapping^for matrices
785 
786 .seealso:  MatAssemblyBegin(), MatAssemblyEnd(), MatSetValues(), MatSetValuesLocal()
787 @*/
788 int MatSetLocalToGlobalMapping(Mat x,ISLocalToGlobalMapping mapping)
789 {
790   int ierr;
791   PetscFunctionBegin;
792   PetscValidHeaderSpecific(x,MAT_COOKIE);
793   PetscValidType(x);
794   MatPreallocated(x);
795   PetscValidHeaderSpecific(mapping,IS_LTOGM_COOKIE);
796   if (x->mapping) {
797     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Mapping already set for matrix");
798   }
799 
800   if (x->ops->setlocaltoglobalmapping) {
801     ierr = (*x->ops->setlocaltoglobalmapping)(x,mapping);CHKERRQ(ierr);
802   } else {
803     x->mapping = mapping;
804     ierr = PetscObjectReference((PetscObject)mapping);CHKERRQ(ierr);
805   }
806   PetscFunctionReturn(0);
807 }
808 
809 #undef __FUNC__
810 #define __FUNC__ "MatSetLocalToGlobalMappingBlock"
811 /*@
812    MatSetLocalToGlobalMappingBlock - Sets a local-to-global numbering for use
813    by the routine MatSetValuesBlockedLocal() to allow users to insert matrix
814    entries using a local (per-processor) numbering.
815 
816    Not Collective
817 
818    Input Parameters:
819 +  x - the matrix
820 -  mapping - mapping created with ISLocalToGlobalMappingCreate() or
821              ISLocalToGlobalMappingCreateIS()
822 
823    Level: intermediate
824 
825    Concepts: matrices^local to global mapping blocked
826    Concepts: local to global mapping^for matrices, blocked
827 
828 .seealso:  MatAssemblyBegin(), MatAssemblyEnd(), MatSetValues(), MatSetValuesBlockedLocal(),
829            MatSetValuesBlocked(), MatSetValuesLocal()
830 @*/
831 int MatSetLocalToGlobalMappingBlock(Mat x,ISLocalToGlobalMapping mapping)
832 {
833   int ierr;
834   PetscFunctionBegin;
835   PetscValidHeaderSpecific(x,MAT_COOKIE);
836   PetscValidType(x);
837   MatPreallocated(x);
838   PetscValidHeaderSpecific(mapping,IS_LTOGM_COOKIE);
839   if (x->bmapping) {
840     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Mapping already set for matrix");
841   }
842 
843   x->bmapping = mapping;
844   ierr = PetscObjectReference((PetscObject)mapping);CHKERRQ(ierr);
845   PetscFunctionReturn(0);
846 }
847 
848 #undef __FUNC__
849 #define __FUNC__ "MatSetValuesLocal"
850 /*@
851    MatSetValuesLocal - Inserts or adds values into certain locations of a matrix,
852    using a local ordering of the nodes.
853 
854    Not Collective
855 
856    Input Parameters:
857 +  x - the matrix
858 .  nrow, irow - number of rows and their local indices
859 .  ncol, icol - number of columns and their local indices
860 .  y -  a logically two-dimensional array of values
861 -  addv - either INSERT_VALUES or ADD_VALUES, where
862    ADD_VALUES adds values to any existing entries, and
863    INSERT_VALUES replaces existing entries with new values
864 
865    Notes:
866    Before calling MatSetValuesLocal(), the user must first set the
867    local-to-global mapping by calling MatSetLocalToGlobalMapping().
868 
869    Calls to MatSetValuesLocal() with the INSERT_VALUES and ADD_VALUES
870    options cannot be mixed without intervening calls to the assembly
871    routines.
872 
873    These values may be cached, so MatAssemblyBegin() and MatAssemblyEnd()
874    MUST be called after all calls to MatSetValuesLocal() have been completed.
875 
876    Level: intermediate
877 
878    Concepts: matrices^putting entries in with local numbering
879 
880 .seealso:  MatAssemblyBegin(), MatAssemblyEnd(), MatSetValues(), MatSetLocalToGlobalMapping(),
881            MatSetValueLocal()
882 @*/
883 int MatSetValuesLocal(Mat mat,int nrow,int *irow,int ncol,int *icol,Scalar *y,InsertMode addv)
884 {
885   int ierr,irowm[2048],icolm[2048];
886 
887   PetscFunctionBegin;
888   PetscValidHeaderSpecific(mat,MAT_COOKIE);
889   PetscValidType(mat);
890   MatPreallocated(mat);
891   PetscValidIntPointer(irow);
892   PetscValidIntPointer(icol);
893   PetscValidScalarPointer(y);
894 
895   if (mat->insertmode == NOT_SET_VALUES) {
896     mat->insertmode = addv;
897   }
898 #if defined(PETSC_USE_BOPT_g)
899   else if (mat->insertmode != addv) {
900     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot mix add values and insert values");
901   }
902   if (!mat->ops->setvalueslocal && (nrow > 2048 || ncol > 2048)) {
903     SETERRQ2(PETSC_ERR_SUP,"Number column/row indices must be <= 2048: are %d %d",nrow,ncol);
904   }
905   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
906 #endif
907 
908   if (mat->assembled) {
909     mat->was_assembled = PETSC_TRUE;
910     mat->assembled     = PETSC_FALSE;
911   }
912   ierr = PetscLogEventBegin(MAT_SetValues,mat,0,0,0);CHKERRQ(ierr);
913   if (!mat->ops->setvalueslocal) {
914     ierr = ISLocalToGlobalMappingApply(mat->mapping,nrow,irow,irowm);CHKERRQ(ierr);
915     ierr = ISLocalToGlobalMappingApply(mat->mapping,ncol,icol,icolm);CHKERRQ(ierr);
916     ierr = (*mat->ops->setvalues)(mat,nrow,irowm,ncol,icolm,y,addv);CHKERRQ(ierr);
917   } else {
918     ierr = (*mat->ops->setvalueslocal)(mat,nrow,irow,ncol,icol,y,addv);CHKERRQ(ierr);
919   }
920   ierr = PetscLogEventEnd(MAT_SetValues,mat,0,0,0);CHKERRQ(ierr);
921   PetscFunctionReturn(0);
922 }
923 
924 #undef __FUNC__
925 #define __FUNC__ "MatSetValuesBlockedLocal"
926 /*@
927    MatSetValuesBlockedLocal - Inserts or adds values into certain locations of a matrix,
928    using a local ordering of the nodes a block at a time.
929 
930    Not Collective
931 
932    Input Parameters:
933 +  x - the matrix
934 .  nrow, irow - number of rows and their local indices
935 .  ncol, icol - number of columns and their local indices
936 .  y -  a logically two-dimensional array of values
937 -  addv - either INSERT_VALUES or ADD_VALUES, where
938    ADD_VALUES adds values to any existing entries, and
939    INSERT_VALUES replaces existing entries with new values
940 
941    Notes:
942    Before calling MatSetValuesBlockedLocal(), the user must first set the
943    local-to-global mapping by calling MatSetLocalToGlobalMappingBlock(),
944    where the mapping MUST be set for matrix blocks, not for matrix elements.
945 
946    Calls to MatSetValuesBlockedLocal() with the INSERT_VALUES and ADD_VALUES
947    options cannot be mixed without intervening calls to the assembly
948    routines.
949 
950    These values may be cached, so MatAssemblyBegin() and MatAssemblyEnd()
951    MUST be called after all calls to MatSetValuesBlockedLocal() have been completed.
952 
953    Level: intermediate
954 
955    Concepts: matrices^putting blocked values in with local numbering
956 
957 .seealso:  MatAssemblyBegin(), MatAssemblyEnd(), MatSetValuesLocal(), MatSetLocalToGlobalMappingBlock(), MatSetValuesBlocked()
958 @*/
959 int MatSetValuesBlockedLocal(Mat mat,int nrow,int *irow,int ncol,int *icol,Scalar *y,InsertMode addv)
960 {
961   int ierr,irowm[2048],icolm[2048];
962 
963   PetscFunctionBegin;
964   PetscValidHeaderSpecific(mat,MAT_COOKIE);
965   PetscValidType(mat);
966   MatPreallocated(mat);
967   PetscValidIntPointer(irow);
968   PetscValidIntPointer(icol);
969   PetscValidScalarPointer(y);
970   if (mat->insertmode == NOT_SET_VALUES) {
971     mat->insertmode = addv;
972   }
973 #if defined(PETSC_USE_BOPT_g)
974   else if (mat->insertmode != addv) {
975     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot mix add values and insert values");
976   }
977   if (!mat->bmapping) {
978     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Local to global never set with MatSetLocalToGlobalMappingBlock()");
979   }
980   if (nrow > 2048 || ncol > 2048) {
981     SETERRQ2(PETSC_ERR_SUP,"Number column/row indices must be <= 2048: are %d %d",nrow,ncol);
982   }
983   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
984 #endif
985 
986   if (mat->assembled) {
987     mat->was_assembled = PETSC_TRUE;
988     mat->assembled     = PETSC_FALSE;
989   }
990   ierr = PetscLogEventBegin(MAT_SetValues,mat,0,0,0);CHKERRQ(ierr);
991   ierr = ISLocalToGlobalMappingApply(mat->bmapping,nrow,irow,irowm);CHKERRQ(ierr);
992   ierr = ISLocalToGlobalMappingApply(mat->bmapping,ncol,icol,icolm);CHKERRQ(ierr);
993   ierr = (*mat->ops->setvaluesblocked)(mat,nrow,irowm,ncol,icolm,y,addv);CHKERRQ(ierr);
994   ierr = PetscLogEventEnd(MAT_SetValues,mat,0,0,0);CHKERRQ(ierr);
995   PetscFunctionReturn(0);
996 }
997 
998 /* --------------------------------------------------------*/
999 #undef __FUNC__
1000 #define __FUNC__ "MatMult"
1001 /*@
1002    MatMult - Computes the matrix-vector product, y = Ax.
1003 
1004    Collective on Mat and Vec
1005 
1006    Input Parameters:
1007 +  mat - the matrix
1008 -  x   - the vector to be multilplied
1009 
1010    Output Parameters:
1011 .  y - the result
1012 
1013    Notes:
1014    The vectors x and y cannot be the same.  I.e., one cannot
1015    call MatMult(A,y,y).
1016 
1017    Level: beginner
1018 
1019    Concepts: matrix-vector product
1020 
1021 .seealso: MatMultTranspose(), MatMultAdd(), MatMultTransposeAdd()
1022 @*/
1023 int MatMult(Mat mat,Vec x,Vec y)
1024 {
1025   int ierr;
1026 
1027   PetscFunctionBegin;
1028   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1029   PetscValidType(mat);
1030   MatPreallocated(mat);
1031   PetscValidHeaderSpecific(x,VEC_COOKIE);
1032   PetscValidHeaderSpecific(y,VEC_COOKIE);
1033 
1034   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
1035   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1036   if (x == y) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"x and y must be different vectors");
1037   if (mat->N != x->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec x: global dim %d %d",mat->N,x->N);
1038   if (mat->M != y->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec y: global dim %d %d",mat->M,y->N);
1039   if (mat->m != y->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec y: local dim %d %d",mat->m,y->n);
1040 
1041   if (mat->nullsp) {
1042     ierr = MatNullSpaceRemove(mat->nullsp,x,&x);CHKERRQ(ierr);
1043   }
1044 
1045   ierr = PetscLogEventBegin(MAT_Mult,mat,x,y,0);CHKERRQ(ierr);
1046   ierr = (*mat->ops->mult)(mat,x,y);CHKERRQ(ierr);
1047   ierr = PetscLogEventEnd(MAT_Mult,mat,x,y,0);CHKERRQ(ierr);
1048 
1049   if (mat->nullsp) {
1050     ierr = MatNullSpaceRemove(mat->nullsp,y,PETSC_NULL);CHKERRQ(ierr);
1051   }
1052   PetscFunctionReturn(0);
1053 }
1054 
1055 #undef __FUNC__
1056 #define __FUNC__ "MatMultTranspose"
1057 /*@
1058    MatMultTranspose - Computes matrix transpose times a vector.
1059 
1060    Collective on Mat and Vec
1061 
1062    Input Parameters:
1063 +  mat - the matrix
1064 -  x   - the vector to be multilplied
1065 
1066    Output Parameters:
1067 .  y - the result
1068 
1069    Notes:
1070    The vectors x and y cannot be the same.  I.e., one cannot
1071    call MatMultTranspose(A,y,y).
1072 
1073    Level: beginner
1074 
1075    Concepts: matrix vector product^transpose
1076 
1077 .seealso: MatMult(), MatMultAdd(), MatMultTransposeAdd()
1078 @*/
1079 int MatMultTranspose(Mat mat,Vec x,Vec y)
1080 {
1081   int ierr;
1082 
1083   PetscFunctionBegin;
1084   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1085   PetscValidType(mat);
1086   MatPreallocated(mat);
1087   PetscValidHeaderSpecific(x,VEC_COOKIE);
1088   PetscValidHeaderSpecific(y,VEC_COOKIE);
1089 
1090   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
1091   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1092   if (x == y) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"x and y must be different vectors");
1093   if (mat->M != x->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec x: global dim %d %d",mat->M,x->N);
1094   if (mat->N != y->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec y: global dim %d %d",mat->N,y->N);
1095 
1096   ierr = PetscLogEventBegin(MAT_MultTranspose,mat,x,y,0);CHKERRQ(ierr);
1097   ierr = (*mat->ops->multtranspose)(mat,x,y);CHKERRQ(ierr);
1098   ierr = PetscLogEventEnd(MAT_MultTranspose,mat,x,y,0);CHKERRQ(ierr);
1099   PetscFunctionReturn(0);
1100 }
1101 
1102 #undef __FUNC__
1103 #define __FUNC__ "MatMultAdd"
1104 /*@
1105     MatMultAdd -  Computes v3 = v2 + A * v1.
1106 
1107     Collective on Mat and Vec
1108 
1109     Input Parameters:
1110 +   mat - the matrix
1111 -   v1, v2 - the vectors
1112 
1113     Output Parameters:
1114 .   v3 - the result
1115 
1116     Notes:
1117     The vectors v1 and v3 cannot be the same.  I.e., one cannot
1118     call MatMultAdd(A,v1,v2,v1).
1119 
1120     Level: beginner
1121 
1122     Concepts: matrix vector product^addition
1123 
1124 .seealso: MatMultTranspose(), MatMult(), MatMultTransposeAdd()
1125 @*/
1126 int MatMultAdd(Mat mat,Vec v1,Vec v2,Vec v3)
1127 {
1128   int ierr;
1129 
1130   PetscFunctionBegin;
1131   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1132   PetscValidType(mat);
1133   MatPreallocated(mat);
1134   PetscValidHeaderSpecific(v1,VEC_COOKIE);
1135   PetscValidHeaderSpecific(v2,VEC_COOKIE);
1136   PetscValidHeaderSpecific(v3,VEC_COOKIE);
1137 
1138   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
1139   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1140   if (mat->N != v1->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec v1: global dim %d %d",mat->N,v1->N);
1141   if (mat->M != v2->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec v2: global dim %d %d",mat->M,v2->N);
1142   if (mat->M != v3->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec v3: global dim %d %d",mat->M,v3->N);
1143   if (mat->m != v3->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec v3: local dim %d %d",mat->m,v3->n);
1144   if (mat->m != v2->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec v2: local dim %d %d",mat->m,v2->n);
1145   if (v1 == v3) SETERRQ(PETSC_ERR_ARG_IDN,"v1 and v3 must be different vectors");
1146 
1147   ierr = PetscLogEventBegin(MAT_MultAdd,mat,v1,v2,v3);CHKERRQ(ierr);
1148   ierr = (*mat->ops->multadd)(mat,v1,v2,v3);CHKERRQ(ierr);
1149   ierr = PetscLogEventEnd(MAT_MultAdd,mat,v1,v2,v3);CHKERRQ(ierr);
1150   PetscFunctionReturn(0);
1151 }
1152 
1153 #undef __FUNC__
1154 #define __FUNC__ "MatMultTransposeAdd"
1155 /*@
1156    MatMultTransposeAdd - Computes v3 = v2 + A' * v1.
1157 
1158    Collective on Mat and Vec
1159 
1160    Input Parameters:
1161 +  mat - the matrix
1162 -  v1, v2 - the vectors
1163 
1164    Output Parameters:
1165 .  v3 - the result
1166 
1167    Notes:
1168    The vectors v1 and v3 cannot be the same.  I.e., one cannot
1169    call MatMultTransposeAdd(A,v1,v2,v1).
1170 
1171    Level: beginner
1172 
1173    Concepts: matrix vector product^transpose and addition
1174 
1175 .seealso: MatMultTranspose(), MatMultAdd(), MatMult()
1176 @*/
1177 int MatMultTransposeAdd(Mat mat,Vec v1,Vec v2,Vec v3)
1178 {
1179   int ierr;
1180 
1181   PetscFunctionBegin;
1182   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1183   PetscValidType(mat);
1184   MatPreallocated(mat);
1185   PetscValidHeaderSpecific(v1,VEC_COOKIE);
1186   PetscValidHeaderSpecific(v2,VEC_COOKIE);
1187   PetscValidHeaderSpecific(v3,VEC_COOKIE);
1188 
1189   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
1190   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1191   if (!mat->ops->multtransposeadd) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name);
1192   if (v1 == v3) SETERRQ(PETSC_ERR_ARG_IDN,"v1 and v3 must be different vectors");
1193   if (mat->M != v1->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec v1: global dim %d %d",mat->M,v1->N);
1194   if (mat->N != v2->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec v2: global dim %d %d",mat->N,v2->N);
1195   if (mat->N != v3->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec v3: global dim %d %d",mat->N,v3->N);
1196 
1197   ierr = PetscLogEventBegin(MAT_MultTransposeAdd,mat,v1,v2,v3);CHKERRQ(ierr);
1198   ierr = (*mat->ops->multtransposeadd)(mat,v1,v2,v3);CHKERRQ(ierr);
1199   ierr = PetscLogEventEnd(MAT_MultTransposeAdd,mat,v1,v2,v3);CHKERRQ(ierr);
1200   PetscFunctionReturn(0);
1201 }
1202 /* ------------------------------------------------------------*/
1203 #undef __FUNC__
1204 #define __FUNC__ "MatGetInfo"
1205 /*@C
1206    MatGetInfo - Returns information about matrix storage (number of
1207    nonzeros, memory, etc.).
1208 
1209    Collective on Mat if MAT_GLOBAL_MAX or MAT_GLOBAL_SUM is used
1210    as the flag
1211 
1212    Input Parameters:
1213 .  mat - the matrix
1214 
1215    Output Parameters:
1216 +  flag - flag indicating the type of parameters to be returned
1217    (MAT_LOCAL - local matrix, MAT_GLOBAL_MAX - maximum over all processors,
1218    MAT_GLOBAL_SUM - sum over all processors)
1219 -  info - matrix information context
1220 
1221    Notes:
1222    The MatInfo context contains a variety of matrix data, including
1223    number of nonzeros allocated and used, number of mallocs during
1224    matrix assembly, etc.  Additional information for factored matrices
1225    is provided (such as the fill ratio, number of mallocs during
1226    factorization, etc.).  Much of this info is printed to STDOUT
1227    when using the runtime options
1228 $       -log_info -mat_view_info
1229 
1230    Example for C/C++ Users:
1231    See the file ${PETSC_DIR}/include/petscmat.h for a complete list of
1232    data within the MatInfo context.  For example,
1233 .vb
1234       MatInfo info;
1235       Mat     A;
1236       double  mal, nz_a, nz_u;
1237 
1238       MatGetInfo(A,MAT_LOCAL,&info);
1239       mal  = info.mallocs;
1240       nz_a = info.nz_allocated;
1241 .ve
1242 
1243    Example for Fortran Users:
1244    Fortran users should declare info as a double precision
1245    array of dimension MAT_INFO_SIZE, and then extract the parameters
1246    of interest.  See the file ${PETSC_DIR}/include/finclude/petscmat.h
1247    a complete list of parameter names.
1248 .vb
1249       double  precision info(MAT_INFO_SIZE)
1250       double  precision mal, nz_a
1251       Mat     A
1252       integer ierr
1253 
1254       call MatGetInfo(A,MAT_LOCAL,info,ierr)
1255       mal = info(MAT_INFO_MALLOCS)
1256       nz_a = info(MAT_INFO_NZ_ALLOCATED)
1257 .ve
1258 
1259     Level: intermediate
1260 
1261     Concepts: matrices^getting information on
1262 
1263 @*/
1264 int MatGetInfo(Mat mat,MatInfoType flag,MatInfo *info)
1265 {
1266   int ierr;
1267 
1268   PetscFunctionBegin;
1269   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1270   PetscValidType(mat);
1271   MatPreallocated(mat);
1272   PetscValidPointer(info);
1273   if (!mat->ops->getinfo) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name);
1274   ierr = (*mat->ops->getinfo)(mat,flag,info);CHKERRQ(ierr);
1275   PetscFunctionReturn(0);
1276 }
1277 
1278 /* ----------------------------------------------------------*/
1279 #undef __FUNC__
1280 #define __FUNC__ "MatILUDTFactor"
1281 /*@C
1282    MatILUDTFactor - Performs a drop tolerance ILU factorization.
1283 
1284    Collective on Mat
1285 
1286    Input Parameters:
1287 +  mat - the matrix
1288 .  info - information about the factorization to be done
1289 .  row - row permutation
1290 -  col - column permutation
1291 
1292    Output Parameters:
1293 .  fact - the factored matrix
1294 
1295    Level: developer
1296 
1297    Notes:
1298    Most users should employ the simplified SLES interface for linear solvers
1299    instead of working directly with matrix algebra routines such as this.
1300    See, e.g., SLESCreate().
1301 
1302    This is currently only supported for the SeqAIJ matrix format using code
1303    from Yousef Saad's SPARSEKIT2  package (translated to C with f2c) and/or
1304    Matlab. SPARSEKIT2 is copyrighted by Yousef Saad with the GNU copyright
1305    and thus can be distributed with PETSc.
1306 
1307     Concepts: matrices^ILUDT factorization
1308 
1309 .seealso: MatLUFactorSymbolic(), MatLUFactorNumeric(), MatCholeskyFactor()
1310 @*/
1311 int MatILUDTFactor(Mat mat,MatILUInfo *info,IS row,IS col,Mat *fact)
1312 {
1313   int ierr;
1314 
1315   PetscFunctionBegin;
1316   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1317   PetscValidType(mat);
1318   MatPreallocated(mat);
1319   PetscValidPointer(fact);
1320   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
1321   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1322   if (!mat->ops->iludtfactor) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name);
1323 
1324   ierr = PetscLogEventBegin(MAT_ILUFactor,mat,row,col,0);CHKERRQ(ierr);
1325   ierr = (*mat->ops->iludtfactor)(mat,info,row,col,fact);CHKERRQ(ierr);
1326   ierr = PetscLogEventEnd(MAT_ILUFactor,mat,row,col,0);CHKERRQ(ierr);
1327 
1328   PetscFunctionReturn(0);
1329 }
1330 
1331 #undef __FUNC__
1332 #define __FUNC__ "MatLUFactor"
1333 /*@
1334    MatLUFactor - Performs in-place LU factorization of matrix.
1335 
1336    Collective on Mat
1337 
1338    Input Parameters:
1339 +  mat - the matrix
1340 .  row - row permutation
1341 .  col - column permutation
1342 -  info - options for factorization, includes
1343 $          fill - expected fill as ratio of original fill.
1344 $          dtcol - pivot tolerance (0 no pivot, 1 full column pivoting)
1345 $                   Run with the option -log_info to determine an optimal value to use
1346 
1347    Notes:
1348    Most users should employ the simplified SLES interface for linear solvers
1349    instead of working directly with matrix algebra routines such as this.
1350    See, e.g., SLESCreate().
1351 
1352    This changes the state of the matrix to a factored matrix; it cannot be used
1353    for example with MatSetValues() unless one first calls MatSetUnfactored().
1354 
1355    Level: developer
1356 
1357    Concepts: matrices^LU factorization
1358 
1359 .seealso: MatLUFactorSymbolic(), MatLUFactorNumeric(), MatCholeskyFactor(),
1360           MatGetOrdering(), MatSetUnfactored()
1361 
1362 @*/
1363 int MatLUFactor(Mat mat,IS row,IS col,MatLUInfo *info)
1364 {
1365   int ierr;
1366 
1367   PetscFunctionBegin;
1368   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1369   PetscValidType(mat);
1370   MatPreallocated(mat);
1371   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
1372   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1373   if (!mat->ops->lufactor) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name);
1374 
1375   ierr = PetscLogEventBegin(MAT_LUFactor,mat,row,col,0);CHKERRQ(ierr);
1376   ierr = (*mat->ops->lufactor)(mat,row,col,info);CHKERRQ(ierr);
1377   ierr = PetscLogEventEnd(MAT_LUFactor,mat,row,col,0);CHKERRQ(ierr);
1378   PetscFunctionReturn(0);
1379 }
1380 
1381 #undef __FUNC__
1382 #define __FUNC__ "MatILUFactor"
1383 /*@
1384    MatILUFactor - Performs in-place ILU factorization of matrix.
1385 
1386    Collective on Mat
1387 
1388    Input Parameters:
1389 +  mat - the matrix
1390 .  row - row permutation
1391 .  col - column permutation
1392 -  info - structure containing
1393 $      levels - number of levels of fill.
1394 $      expected fill - as ratio of original fill.
1395 $      1 or 0 - indicating force fill on diagonal (improves robustness for matrices
1396                 missing diagonal entries)
1397 
1398    Notes:
1399    Probably really in-place only when level of fill is zero, otherwise allocates
1400    new space to store factored matrix and deletes previous memory.
1401 
1402    Most users should employ the simplified SLES interface for linear solvers
1403    instead of working directly with matrix algebra routines such as this.
1404    See, e.g., SLESCreate().
1405 
1406    Level: developer
1407 
1408    Concepts: matrices^ILU factorization
1409 
1410 .seealso: MatILUFactorSymbolic(), MatLUFactorNumeric(), MatCholeskyFactor()
1411 @*/
1412 int MatILUFactor(Mat mat,IS row,IS col,MatILUInfo *info)
1413 {
1414   int ierr;
1415 
1416   PetscFunctionBegin;
1417   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1418   PetscValidType(mat);
1419   MatPreallocated(mat);
1420   if (mat->M != mat->N) SETERRQ(PETSC_ERR_ARG_WRONG,"matrix must be square");
1421   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
1422   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1423   if (!mat->ops->ilufactor) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name);
1424 
1425   ierr = PetscLogEventBegin(MAT_ILUFactor,mat,row,col,0);CHKERRQ(ierr);
1426   ierr = (*mat->ops->ilufactor)(mat,row,col,info);CHKERRQ(ierr);
1427   ierr = PetscLogEventEnd(MAT_ILUFactor,mat,row,col,0);CHKERRQ(ierr);
1428   PetscFunctionReturn(0);
1429 }
1430 
1431 #undef __FUNC__
1432 #define __FUNC__ "MatLUFactorSymbolic"
1433 /*@
1434    MatLUFactorSymbolic - Performs symbolic LU factorization of matrix.
1435    Call this routine before calling MatLUFactorNumeric().
1436 
1437    Collective on Mat
1438 
1439    Input Parameters:
1440 +  mat - the matrix
1441 .  row, col - row and column permutations
1442 -  info - options for factorization, includes
1443 $          fill - expected fill as ratio of original fill.
1444 $          dtcol - pivot tolerance (0 no pivot, 1 full column pivoting)
1445 $                   Run with the option -log_info to determine an optimal value to use
1446 
1447    Output Parameter:
1448 .  fact - new matrix that has been symbolically factored
1449 
1450    Notes:
1451    See the users manual for additional information about
1452    choosing the fill factor for better efficiency.
1453 
1454    Most users should employ the simplified SLES interface for linear solvers
1455    instead of working directly with matrix algebra routines such as this.
1456    See, e.g., SLESCreate().
1457 
1458    Level: developer
1459 
1460    Concepts: matrices^LU symbolic factorization
1461 
1462 .seealso: MatLUFactor(), MatLUFactorNumeric(), MatCholeskyFactor()
1463 @*/
1464 int MatLUFactorSymbolic(Mat mat,IS row,IS col,MatLUInfo *info,Mat *fact)
1465 {
1466   int ierr;
1467 
1468   PetscFunctionBegin;
1469   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1470   PetscValidType(mat);
1471   MatPreallocated(mat);
1472   PetscValidPointer(fact);
1473   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
1474   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1475   if (!mat->ops->lufactorsymbolic) SETERRQ1(PETSC_ERR_SUP,"Matrix type %s  symbolic LU",mat->type_name);
1476 
1477   ierr = PetscLogEventBegin(MAT_LUFactorSymbolic,mat,row,col,0);CHKERRQ(ierr);
1478   ierr = (*mat->ops->lufactorsymbolic)(mat,row,col,info,fact);CHKERRQ(ierr);
1479   ierr = PetscLogEventEnd(MAT_LUFactorSymbolic,mat,row,col,0);CHKERRQ(ierr);
1480   PetscFunctionReturn(0);
1481 }
1482 
1483 #undef __FUNC__
1484 #define __FUNC__ "MatLUFactorNumeric"
1485 /*@
1486    MatLUFactorNumeric - Performs numeric LU factorization of a matrix.
1487    Call this routine after first calling MatLUFactorSymbolic().
1488 
1489    Collective on Mat
1490 
1491    Input Parameters:
1492 +  mat - the matrix
1493 -  fact - the matrix generated for the factor, from MatLUFactorSymbolic()
1494 
1495    Notes:
1496    See MatLUFactor() for in-place factorization.  See
1497    MatCholeskyFactorNumeric() for the symmetric, positive definite case.
1498 
1499    Most users should employ the simplified SLES interface for linear solvers
1500    instead of working directly with matrix algebra routines such as this.
1501    See, e.g., SLESCreate().
1502 
1503    Level: developer
1504 
1505    Concepts: matrices^LU numeric factorization
1506 
1507 .seealso: MatLUFactorSymbolic(), MatLUFactor(), MatCholeskyFactor()
1508 @*/
1509 int MatLUFactorNumeric(Mat mat,Mat *fact)
1510 {
1511   int        ierr;
1512   PetscTruth flg;
1513 
1514   PetscFunctionBegin;
1515   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1516   PetscValidType(mat);
1517   MatPreallocated(mat);
1518   PetscValidPointer(fact);
1519   PetscValidHeaderSpecific(*fact,MAT_COOKIE);
1520   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
1521   if (mat->M != (*fact)->M || mat->N != (*fact)->N) {
1522     SETERRQ4(PETSC_ERR_ARG_SIZ,"Mat mat,Mat *fact: global dimensions are different %d should = %d %d should = %d",
1523             mat->M,(*fact)->M,mat->N,(*fact)->N);
1524   }
1525   if (!(*fact)->ops->lufactornumeric) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name);
1526 
1527   ierr = PetscLogEventBegin(MAT_LUFactorNumeric,mat,*fact,0,0);CHKERRQ(ierr);
1528   ierr = (*(*fact)->ops->lufactornumeric)(mat,fact);CHKERRQ(ierr);
1529   ierr = PetscLogEventEnd(MAT_LUFactorNumeric,mat,*fact,0,0);CHKERRQ(ierr);
1530   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_view_draw",&flg);CHKERRQ(ierr);
1531   if (flg) {
1532     ierr = PetscOptionsHasName(PETSC_NULL,"-mat_view_contour",&flg);CHKERRQ(ierr);
1533     if (flg) {
1534       ierr = PetscViewerPushFormat(PETSC_VIEWER_DRAW_(mat->comm),PETSC_VIEWER_DRAW_CONTOUR);CHKERRQ(ierr);
1535     }
1536     ierr = MatView(*fact,PETSC_VIEWER_DRAW_(mat->comm));CHKERRQ(ierr);
1537     ierr = PetscViewerFlush(PETSC_VIEWER_DRAW_(mat->comm));CHKERRQ(ierr);
1538     if (flg) {
1539       ierr = PetscViewerPopFormat(PETSC_VIEWER_DRAW_(mat->comm));CHKERRQ(ierr);
1540     }
1541   }
1542   PetscFunctionReturn(0);
1543 }
1544 
1545 #undef __FUNC__
1546 #define __FUNC__ "MatCholeskyFactor"
1547 /*@
1548    MatCholeskyFactor - Performs in-place Cholesky factorization of a
1549    symmetric matrix.
1550 
1551    Collective on Mat
1552 
1553    Input Parameters:
1554 +  mat - the matrix
1555 .  perm - row and column permutations
1556 -  f - expected fill as ratio of original fill
1557 
1558    Notes:
1559    See MatLUFactor() for the nonsymmetric case.  See also
1560    MatCholeskyFactorSymbolic(), and MatCholeskyFactorNumeric().
1561 
1562    Most users should employ the simplified SLES interface for linear solvers
1563    instead of working directly with matrix algebra routines such as this.
1564    See, e.g., SLESCreate().
1565 
1566    Level: developer
1567 
1568    Concepts: matrices^Cholesky factorization
1569 
1570 .seealso: MatLUFactor(), MatCholeskyFactorSymbolic(), MatCholeskyFactorNumeric()
1571           MatGetOrdering()
1572 
1573 @*/
1574 int MatCholeskyFactor(Mat mat,IS perm,PetscReal f)
1575 {
1576   int ierr;
1577 
1578   PetscFunctionBegin;
1579   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1580   PetscValidType(mat);
1581   MatPreallocated(mat);
1582   if (mat->M != mat->N) SETERRQ(PETSC_ERR_ARG_WRONG,"Matrix must be square");
1583   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
1584   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1585   if (!mat->ops->choleskyfactor) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name);
1586 
1587   ierr = PetscLogEventBegin(MAT_CholeskyFactor,mat,perm,0,0);CHKERRQ(ierr);
1588   ierr = (*mat->ops->choleskyfactor)(mat,perm,f);CHKERRQ(ierr);
1589   ierr = PetscLogEventEnd(MAT_CholeskyFactor,mat,perm,0,0);CHKERRQ(ierr);
1590   PetscFunctionReturn(0);
1591 }
1592 
1593 #undef __FUNC__
1594 #define __FUNC__ "MatCholeskyFactorSymbolic"
1595 /*@
1596    MatCholeskyFactorSymbolic - Performs symbolic Cholesky factorization
1597    of a symmetric matrix.
1598 
1599    Collective on Mat
1600 
1601    Input Parameters:
1602 +  mat - the matrix
1603 .  perm - row and column permutations
1604 -  f - expected fill as ratio of original
1605 
1606    Output Parameter:
1607 .  fact - the factored matrix
1608 
1609    Notes:
1610    See MatLUFactorSymbolic() for the nonsymmetric case.  See also
1611    MatCholeskyFactor() and MatCholeskyFactorNumeric().
1612 
1613    Most users should employ the simplified SLES interface for linear solvers
1614    instead of working directly with matrix algebra routines such as this.
1615    See, e.g., SLESCreate().
1616 
1617    Level: developer
1618 
1619    Concepts: matrices^Cholesky symbolic factorization
1620 
1621 .seealso: MatLUFactorSymbolic(), MatCholeskyFactor(), MatCholeskyFactorNumeric()
1622           MatGetOrdering()
1623 
1624 @*/
1625 int MatCholeskyFactorSymbolic(Mat mat,IS perm,PetscReal f,Mat *fact)
1626 {
1627   int ierr;
1628 
1629   PetscFunctionBegin;
1630   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1631   PetscValidType(mat);
1632   MatPreallocated(mat);
1633   PetscValidPointer(fact);
1634   if (mat->M != mat->N) SETERRQ(PETSC_ERR_ARG_WRONG,"Matrix must be square");
1635   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
1636   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1637   if (!mat->ops->choleskyfactorsymbolic) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name);
1638 
1639   ierr = PetscLogEventBegin(MAT_CholeskyFactorSymbolic,mat,perm,0,0);CHKERRQ(ierr);
1640   ierr = (*mat->ops->choleskyfactorsymbolic)(mat,perm,f,fact);CHKERRQ(ierr);
1641   ierr = PetscLogEventEnd(MAT_CholeskyFactorSymbolic,mat,perm,0,0);CHKERRQ(ierr);
1642   PetscFunctionReturn(0);
1643 }
1644 
1645 #undef __FUNC__
1646 #define __FUNC__ "MatCholeskyFactorNumeric"
1647 /*@
1648    MatCholeskyFactorNumeric - Performs numeric Cholesky factorization
1649    of a symmetric matrix. Call this routine after first calling
1650    MatCholeskyFactorSymbolic().
1651 
1652    Collective on Mat
1653 
1654    Input Parameter:
1655 .  mat - the initial matrix
1656 
1657    Output Parameter:
1658 .  fact - the factored matrix
1659 
1660    Notes:
1661    Most users should employ the simplified SLES interface for linear solvers
1662    instead of working directly with matrix algebra routines such as this.
1663    See, e.g., SLESCreate().
1664 
1665    Level: developer
1666 
1667    Concepts: matrices^Cholesky numeric factorization
1668 
1669 .seealso: MatCholeskyFactorSymbolic(), MatCholeskyFactor(), MatLUFactorNumeric()
1670 @*/
1671 int MatCholeskyFactorNumeric(Mat mat,Mat *fact)
1672 {
1673   int ierr;
1674 
1675   PetscFunctionBegin;
1676   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1677   PetscValidType(mat);
1678   MatPreallocated(mat);
1679   PetscValidPointer(fact);
1680   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
1681   if (!(*fact)->ops->choleskyfactornumeric) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name);
1682   if (mat->M != (*fact)->M || mat->N != (*fact)->N) {
1683     SETERRQ4(PETSC_ERR_ARG_SIZ,"Mat mat,Mat *fact: global dim %d should = %d %d should = %d",
1684             mat->M,(*fact)->M,mat->N,(*fact)->N);
1685   }
1686 
1687   ierr = PetscLogEventBegin(MAT_CholeskyFactorNumeric,mat,*fact,0,0);CHKERRQ(ierr);
1688   ierr = (*(*fact)->ops->choleskyfactornumeric)(mat,fact);CHKERRQ(ierr);
1689   ierr = PetscLogEventEnd(MAT_CholeskyFactorNumeric,mat,*fact,0,0);CHKERRQ(ierr);
1690   PetscFunctionReturn(0);
1691 }
1692 
1693 /* ----------------------------------------------------------------*/
1694 #undef __FUNC__
1695 #define __FUNC__ "MatSolve"
1696 /*@
1697    MatSolve - Solves A x = b, given a factored matrix.
1698 
1699    Collective on Mat and Vec
1700 
1701    Input Parameters:
1702 +  mat - the factored matrix
1703 -  b - the right-hand-side vector
1704 
1705    Output Parameter:
1706 .  x - the result vector
1707 
1708    Notes:
1709    The vectors b and x cannot be the same.  I.e., one cannot
1710    call MatSolve(A,x,x).
1711 
1712    Notes:
1713    Most users should employ the simplified SLES interface for linear solvers
1714    instead of working directly with matrix algebra routines such as this.
1715    See, e.g., SLESCreate().
1716 
1717    Level: developer
1718 
1719    Concepts: matrices^triangular solves
1720 
1721 .seealso: MatSolveAdd(), MatSolveTranspose(), MatSolveTransposeAdd()
1722 @*/
1723 int MatSolve(Mat mat,Vec b,Vec x)
1724 {
1725   int ierr;
1726 
1727   PetscFunctionBegin;
1728   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1729   PetscValidType(mat);
1730   MatPreallocated(mat);
1731   PetscValidHeaderSpecific(b,VEC_COOKIE);
1732   PetscValidHeaderSpecific(x,VEC_COOKIE);
1733   PetscCheckSameComm(mat,b);
1734   PetscCheckSameComm(mat,x);
1735   if (x == b) SETERRQ(PETSC_ERR_ARG_IDN,"x and b must be different vectors");
1736   if (!mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Unfactored matrix");
1737   if (mat->N != x->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec x: global dim %d %d",mat->N,x->N);
1738   if (mat->M != b->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec b: global dim %d %d",mat->M,b->N);
1739   if (mat->m != b->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec b: local dim %d %d",mat->m,b->n);
1740   if (mat->M == 0 && mat->N == 0) PetscFunctionReturn(0);
1741 
1742   if (!mat->ops->solve) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name);
1743   ierr = PetscLogEventBegin(MAT_Solve,mat,b,x,0);CHKERRQ(ierr);
1744   ierr = (*mat->ops->solve)(mat,b,x);CHKERRQ(ierr);
1745   ierr = PetscLogEventEnd(MAT_Solve,mat,b,x,0);CHKERRQ(ierr);
1746   PetscFunctionReturn(0);
1747 }
1748 
1749 #undef __FUNC__
1750 #define __FUNC__ "MatForwardSolve"
1751 /* @
1752    MatForwardSolve - Solves L x = b, given a factored matrix, A = LU.
1753 
1754    Collective on Mat and Vec
1755 
1756    Input Parameters:
1757 +  mat - the factored matrix
1758 -  b - the right-hand-side vector
1759 
1760    Output Parameter:
1761 .  x - the result vector
1762 
1763    Notes:
1764    MatSolve() should be used for most applications, as it performs
1765    a forward solve followed by a backward solve.
1766 
1767    The vectors b and x cannot be the same.  I.e., one cannot
1768    call MatForwardSolve(A,x,x).
1769 
1770    Most users should employ the simplified SLES interface for linear solvers
1771    instead of working directly with matrix algebra routines such as this.
1772    See, e.g., SLESCreate().
1773 
1774    Level: developer
1775 
1776    Concepts: matrices^forward solves
1777 
1778 .seealso: MatSolve(), MatBackwardSolve()
1779 @ */
1780 int MatForwardSolve(Mat mat,Vec b,Vec x)
1781 {
1782   int ierr;
1783 
1784   PetscFunctionBegin;
1785   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1786   PetscValidType(mat);
1787   MatPreallocated(mat);
1788   PetscValidHeaderSpecific(b,VEC_COOKIE);
1789   PetscValidHeaderSpecific(x,VEC_COOKIE);
1790   PetscCheckSameComm(mat,b);
1791   PetscCheckSameComm(mat,x);
1792   if (x == b) SETERRQ(PETSC_ERR_ARG_IDN,"x and b must be different vectors");
1793   if (!mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Unfactored matrix");
1794   if (!mat->ops->forwardsolve) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name);
1795   if (mat->N != x->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec x: global dim %d %d",mat->N,x->N);
1796   if (mat->M != b->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec b: global dim %d %d",mat->M,b->N);
1797   if (mat->m != b->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec b: local dim %d %d",mat->m,b->n);
1798 
1799   ierr = PetscLogEventBegin(MAT_ForwardSolve,mat,b,x,0);CHKERRQ(ierr);
1800   ierr = (*mat->ops->forwardsolve)(mat,b,x);CHKERRQ(ierr);
1801   ierr = PetscLogEventEnd(MAT_ForwardSolve,mat,b,x,0);CHKERRQ(ierr);
1802   PetscFunctionReturn(0);
1803 }
1804 
1805 #undef __FUNC__
1806 #define __FUNC__ "MatBackwardSolve"
1807 /* @
1808    MatBackwardSolve - Solves U x = b, given a factored matrix, A = LU.
1809 
1810    Collective on Mat and Vec
1811 
1812    Input Parameters:
1813 +  mat - the factored matrix
1814 -  b - the right-hand-side vector
1815 
1816    Output Parameter:
1817 .  x - the result vector
1818 
1819    Notes:
1820    MatSolve() should be used for most applications, as it performs
1821    a forward solve followed by a backward solve.
1822 
1823    The vectors b and x cannot be the same.  I.e., one cannot
1824    call MatBackwardSolve(A,x,x).
1825 
1826    Most users should employ the simplified SLES interface for linear solvers
1827    instead of working directly with matrix algebra routines such as this.
1828    See, e.g., SLESCreate().
1829 
1830    Level: developer
1831 
1832    Concepts: matrices^backward solves
1833 
1834 .seealso: MatSolve(), MatForwardSolve()
1835 @ */
1836 int MatBackwardSolve(Mat mat,Vec b,Vec x)
1837 {
1838   int ierr;
1839 
1840   PetscFunctionBegin;
1841   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1842   PetscValidType(mat);
1843   MatPreallocated(mat);
1844   PetscValidHeaderSpecific(b,VEC_COOKIE);
1845   PetscValidHeaderSpecific(x,VEC_COOKIE);
1846   PetscCheckSameComm(mat,b);
1847   PetscCheckSameComm(mat,x);
1848   if (x == b) SETERRQ(PETSC_ERR_ARG_IDN,"x and b must be different vectors");
1849   if (!mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Unfactored matrix");
1850   if (!mat->ops->backwardsolve) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name);
1851   if (mat->N != x->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec x: global dim %d %d",mat->N,x->N);
1852   if (mat->M != b->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec b: global dim %d %d",mat->M,b->N);
1853   if (mat->m != b->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec b: local dim %d %d",mat->m,b->n);
1854 
1855   ierr = PetscLogEventBegin(MAT_BackwardSolve,mat,b,x,0);CHKERRQ(ierr);
1856   ierr = (*mat->ops->backwardsolve)(mat,b,x);CHKERRQ(ierr);
1857   ierr = PetscLogEventEnd(MAT_BackwardSolve,mat,b,x,0);CHKERRQ(ierr);
1858   PetscFunctionReturn(0);
1859 }
1860 
1861 #undef __FUNC__
1862 #define __FUNC__ "MatSolveAdd"
1863 /*@
1864    MatSolveAdd - Computes x = y + inv(A)*b, given a factored matrix.
1865 
1866    Collective on Mat and Vec
1867 
1868    Input Parameters:
1869 +  mat - the factored matrix
1870 .  b - the right-hand-side vector
1871 -  y - the vector to be added to
1872 
1873    Output Parameter:
1874 .  x - the result vector
1875 
1876    Notes:
1877    The vectors b and x cannot be the same.  I.e., one cannot
1878    call MatSolveAdd(A,x,y,x).
1879 
1880    Most users should employ the simplified SLES interface for linear solvers
1881    instead of working directly with matrix algebra routines such as this.
1882    See, e.g., SLESCreate().
1883 
1884    Level: developer
1885 
1886    Concepts: matrices^triangular solves
1887 
1888 .seealso: MatSolve(), MatSolveTranspose(), MatSolveTransposeAdd()
1889 @*/
1890 int MatSolveAdd(Mat mat,Vec b,Vec y,Vec x)
1891 {
1892   Scalar one = 1.0;
1893   Vec    tmp;
1894   int    ierr;
1895 
1896   PetscFunctionBegin;
1897   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1898   PetscValidType(mat);
1899   MatPreallocated(mat);
1900   PetscValidHeaderSpecific(y,VEC_COOKIE);
1901   PetscValidHeaderSpecific(b,VEC_COOKIE);
1902   PetscValidHeaderSpecific(x,VEC_COOKIE);
1903   PetscCheckSameComm(mat,b);
1904   PetscCheckSameComm(mat,y);
1905   PetscCheckSameComm(mat,x);
1906   if (x == b) SETERRQ(PETSC_ERR_ARG_IDN,"x and b must be different vectors");
1907   if (!mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Unfactored matrix");
1908   if (mat->N != x->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec x: global dim %d %d",mat->N,x->N);
1909   if (mat->M != b->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec b: global dim %d %d",mat->M,b->N);
1910   if (mat->M != y->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec y: global dim %d %d",mat->M,y->N);
1911   if (mat->m != b->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec b: local dim %d %d",mat->m,b->n);
1912   if (x->n != y->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Vec x,Vec y: local dim %d %d",x->n,y->n);
1913 
1914   ierr = PetscLogEventBegin(MAT_SolveAdd,mat,b,x,y);CHKERRQ(ierr);
1915   if (mat->ops->solveadd)  {
1916     ierr = (*mat->ops->solveadd)(mat,b,y,x);CHKERRQ(ierr);
1917   } else {
1918     /* do the solve then the add manually */
1919     if (x != y) {
1920       ierr = MatSolve(mat,b,x);CHKERRQ(ierr);
1921       ierr = VecAXPY(&one,y,x);CHKERRQ(ierr);
1922     } else {
1923       ierr = VecDuplicate(x,&tmp);CHKERRQ(ierr);
1924       PetscLogObjectParent(mat,tmp);
1925       ierr = VecCopy(x,tmp);CHKERRQ(ierr);
1926       ierr = MatSolve(mat,b,x);CHKERRQ(ierr);
1927       ierr = VecAXPY(&one,tmp,x);CHKERRQ(ierr);
1928       ierr = VecDestroy(tmp);CHKERRQ(ierr);
1929     }
1930   }
1931   ierr = PetscLogEventEnd(MAT_SolveAdd,mat,b,x,y);CHKERRQ(ierr);
1932   PetscFunctionReturn(0);
1933 }
1934 
1935 #undef __FUNC__
1936 #define __FUNC__ "MatSolveTranspose"
1937 /*@
1938    MatSolveTranspose - Solves A' x = b, given a factored matrix.
1939 
1940    Collective on Mat and Vec
1941 
1942    Input Parameters:
1943 +  mat - the factored matrix
1944 -  b - the right-hand-side vector
1945 
1946    Output Parameter:
1947 .  x - the result vector
1948 
1949    Notes:
1950    The vectors b and x cannot be the same.  I.e., one cannot
1951    call MatSolveTranspose(A,x,x).
1952 
1953    Most users should employ the simplified SLES interface for linear solvers
1954    instead of working directly with matrix algebra routines such as this.
1955    See, e.g., SLESCreate().
1956 
1957    Level: developer
1958 
1959    Concepts: matrices^triangular solves
1960 
1961 .seealso: MatSolve(), MatSolveAdd(), MatSolveTransposeAdd()
1962 @*/
1963 int MatSolveTranspose(Mat mat,Vec b,Vec x)
1964 {
1965   int ierr;
1966 
1967   PetscFunctionBegin;
1968   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1969   PetscValidType(mat);
1970   MatPreallocated(mat);
1971   PetscValidHeaderSpecific(b,VEC_COOKIE);
1972   PetscValidHeaderSpecific(x,VEC_COOKIE);
1973   PetscCheckSameComm(mat,b);
1974   PetscCheckSameComm(mat,x);
1975   if (!mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Unfactored matrix");
1976   if (x == b) SETERRQ(PETSC_ERR_ARG_IDN,"x and b must be different vectors");
1977   if (!mat->ops->solvetranspose) SETERRQ1(PETSC_ERR_SUP,"Matrix type %s",mat->type_name);
1978   if (mat->M != x->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec x: global dim %d %d",mat->M,x->N);
1979   if (mat->N != b->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec b: global dim %d %d",mat->N,b->N);
1980 
1981   ierr = PetscLogEventBegin(MAT_SolveTranspose,mat,b,x,0);CHKERRQ(ierr);
1982   ierr = (*mat->ops->solvetranspose)(mat,b,x);CHKERRQ(ierr);
1983   ierr = PetscLogEventEnd(MAT_SolveTranspose,mat,b,x,0);CHKERRQ(ierr);
1984   PetscFunctionReturn(0);
1985 }
1986 
1987 #undef __FUNC__
1988 #define __FUNC__ "MatSolveTransposeAdd"
1989 /*@
1990    MatSolveTransposeAdd - Computes x = y + inv(Transpose(A)) b, given a
1991                       factored matrix.
1992 
1993    Collective on Mat and Vec
1994 
1995    Input Parameters:
1996 +  mat - the factored matrix
1997 .  b - the right-hand-side vector
1998 -  y - the vector to be added to
1999 
2000    Output Parameter:
2001 .  x - the result vector
2002 
2003    Notes:
2004    The vectors b and x cannot be the same.  I.e., one cannot
2005    call MatSolveTransposeAdd(A,x,y,x).
2006 
2007    Most users should employ the simplified SLES interface for linear solvers
2008    instead of working directly with matrix algebra routines such as this.
2009    See, e.g., SLESCreate().
2010 
2011    Level: developer
2012 
2013    Concepts: matrices^triangular solves
2014 
2015 .seealso: MatSolve(), MatSolveAdd(), MatSolveTranspose()
2016 @*/
2017 int MatSolveTransposeAdd(Mat mat,Vec b,Vec y,Vec x)
2018 {
2019   Scalar one = 1.0;
2020   int    ierr;
2021   Vec    tmp;
2022 
2023   PetscFunctionBegin;
2024   PetscValidHeaderSpecific(mat,MAT_COOKIE);
2025   PetscValidType(mat);
2026   MatPreallocated(mat);
2027   PetscValidHeaderSpecific(y,VEC_COOKIE);
2028   PetscValidHeaderSpecific(b,VEC_COOKIE);
2029   PetscValidHeaderSpecific(x,VEC_COOKIE);
2030   PetscCheckSameComm(mat,b);
2031   PetscCheckSameComm(mat,y);
2032   PetscCheckSameComm(mat,x);
2033   if (x == b) SETERRQ(PETSC_ERR_ARG_IDN,"x and b must be different vectors");
2034   if (!mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Unfactored matrix");
2035   if (mat->M != x->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec x: global dim %d %d",mat->M,x->N);
2036   if (mat->N != b->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec b: global dim %d %d",mat->N,b->N);
2037   if (mat->N != y->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec y: global dim %d %d",mat->N,y->N);
2038   if (x->n != y->n)   SETERRQ2(PETSC_ERR_ARG_SIZ,"Vec x,Vec y: local dim %d %d",x->n,y->n);
2039 
2040   ierr = PetscLogEventBegin(MAT_SolveTransposeAdd,mat,b,x,y);CHKERRQ(ierr);
2041   if (mat->ops->solvetransposeadd) {
2042     ierr = (*mat->ops->solvetransposeadd)(mat,b,y,x);CHKERRQ(ierr);
2043   } else {
2044     /* do the solve then the add manually */
2045     if (x != y) {
2046       ierr = MatSolveTranspose(mat,b,x);CHKERRQ(ierr);
2047       ierr = VecAXPY(&one,y,x);CHKERRQ(ierr);
2048     } else {
2049       ierr = VecDuplicate(x,&tmp);CHKERRQ(ierr);
2050       PetscLogObjectParent(mat,tmp);
2051       ierr = VecCopy(x,tmp);CHKERRQ(ierr);
2052       ierr = MatSolveTranspose(mat,b,x);CHKERRQ(ierr);
2053       ierr = VecAXPY(&one,tmp,x);CHKERRQ(ierr);
2054       ierr = VecDestroy(tmp);CHKERRQ(ierr);
2055     }
2056   }
2057   ierr = PetscLogEventEnd(MAT_SolveTransposeAdd,mat,b,x,y);CHKERRQ(ierr);
2058   PetscFunctionReturn(0);
2059 }
2060 /* ----------------------------------------------------------------*/
2061 
2062 #undef __FUNC__
2063 #define __FUNC__ "MatRelax"
2064 /*@
2065    MatRelax - Computes one relaxation sweep.
2066 
2067    Collective on Mat and Vec
2068 
2069    Input Parameters:
2070 +  mat - the matrix
2071 .  b - the right hand side
2072 .  omega - the relaxation factor
2073 .  flag - flag indicating the type of SOR (see below)
2074 .  shift -  diagonal shift
2075 -  its - the number of iterations
2076 
2077    Output Parameters:
2078 .  x - the solution (can contain an initial guess)
2079 
2080    SOR Flags:
2081 .     SOR_FORWARD_SWEEP - forward SOR
2082 .     SOR_BACKWARD_SWEEP - backward SOR
2083 .     SOR_SYMMETRIC_SWEEP - SSOR (symmetric SOR)
2084 .     SOR_LOCAL_FORWARD_SWEEP - local forward SOR
2085 .     SOR_LOCAL_BACKWARD_SWEEP - local forward SOR
2086 .     SOR_LOCAL_SYMMETRIC_SWEEP - local SSOR
2087 .     SOR_APPLY_UPPER, SOR_APPLY_LOWER - applies
2088          upper/lower triangular part of matrix to
2089          vector (with omega)
2090 .     SOR_ZERO_INITIAL_GUESS - zero initial guess
2091 
2092    Notes:
2093    SOR_LOCAL_FORWARD_SWEEP, SOR_LOCAL_BACKWARD_SWEEP, and
2094    SOR_LOCAL_SYMMETRIC_SWEEP perform seperate independent smoothings
2095    on each processor.
2096 
2097    Application programmers will not generally use MatRelax() directly,
2098    but instead will employ the SLES/PC interface.
2099 
2100    Notes for Advanced Users:
2101    The flags are implemented as bitwise inclusive or operations.
2102    For example, use (SOR_ZERO_INITIAL_GUESS | SOR_SYMMETRIC_SWEEP)
2103    to specify a zero initial guess for SSOR.
2104 
2105    Most users should employ the simplified SLES interface for linear solvers
2106    instead of working directly with matrix algebra routines such as this.
2107    See, e.g., SLESCreate().
2108 
2109    Level: developer
2110 
2111    Concepts: matrices^relaxation
2112    Concepts: matrices^SOR
2113    Concepts: matrices^Gauss-Seidel
2114 
2115 @*/
2116 int MatRelax(Mat mat,Vec b,PetscReal omega,MatSORType flag,PetscReal shift,int its,Vec x)
2117 {
2118   int ierr;
2119 
2120   PetscFunctionBegin;
2121   PetscValidHeaderSpecific(mat,MAT_COOKIE);
2122   PetscValidType(mat);
2123   MatPreallocated(mat);
2124   PetscValidHeaderSpecific(b,VEC_COOKIE);
2125   PetscValidHeaderSpecific(x,VEC_COOKIE);
2126   PetscCheckSameComm(mat,b);
2127   PetscCheckSameComm(mat,x);
2128   if (!mat->ops->relax) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name);
2129   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2130   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2131   if (mat->N != x->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec x: global dim %d %d",mat->N,x->N);
2132   if (mat->M != b->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec b: global dim %d %d",mat->M,b->N);
2133   if (mat->m != b->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Mat mat,Vec b: local dim %d %d",mat->m,b->n);
2134 
2135   ierr = PetscLogEventBegin(MAT_Relax,mat,b,x,0);CHKERRQ(ierr);
2136   ierr =(*mat->ops->relax)(mat,b,omega,flag,shift,its,x);CHKERRQ(ierr);
2137   ierr = PetscLogEventEnd(MAT_Relax,mat,b,x,0);CHKERRQ(ierr);
2138   PetscFunctionReturn(0);
2139 }
2140 
2141 #undef __FUNC__
2142 #define __FUNC__ "MatCopy_Basic"
2143 /*
2144       Default matrix copy routine.
2145 */
2146 int MatCopy_Basic(Mat A,Mat B,MatStructure str)
2147 {
2148   int    ierr,i,rstart,rend,nz,*cwork;
2149   Scalar *vwork;
2150 
2151   PetscFunctionBegin;
2152   ierr = MatZeroEntries(B);CHKERRQ(ierr);
2153   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
2154   for (i=rstart; i<rend; i++) {
2155     ierr = MatGetRow(A,i,&nz,&cwork,&vwork);CHKERRQ(ierr);
2156     ierr = MatSetValues(B,1,&i,nz,cwork,vwork,INSERT_VALUES);CHKERRQ(ierr);
2157     ierr = MatRestoreRow(A,i,&nz,&cwork,&vwork);CHKERRQ(ierr);
2158   }
2159   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2160   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2161   PetscFunctionReturn(0);
2162 }
2163 
2164 #undef __FUNC__
2165 #define __FUNC__ "MatCopy"
2166 /*@C
2167    MatCopy - Copys a matrix to another matrix.
2168 
2169    Collective on Mat
2170 
2171    Input Parameters:
2172 +  A - the matrix
2173 -  str - SAME_NONZERO_PATTERN or DIFFERENT_NONZERO_PATTERN
2174 
2175    Output Parameter:
2176 .  B - where the copy is put
2177 
2178    Notes:
2179    If you use SAME_NONZERO_PATTERN then the zero matrices had better have the
2180    same nonzero pattern or the routine will crash.
2181 
2182    MatCopy() copies the matrix entries of a matrix to another existing
2183    matrix (after first zeroing the second matrix).  A related routine is
2184    MatConvert(), which first creates a new matrix and then copies the data.
2185 
2186    Level: intermediate
2187 
2188    Concepts: matrices^copying
2189 
2190 .seealso: MatConvert()
2191 @*/
2192 int MatCopy(Mat A,Mat B,MatStructure str)
2193 {
2194   int ierr;
2195 
2196   PetscFunctionBegin;
2197   PetscValidHeaderSpecific(A,MAT_COOKIE);
2198   PetscValidHeaderSpecific(B,MAT_COOKIE);
2199   PetscValidType(A);
2200   MatPreallocated(A);
2201   PetscValidType(B);
2202   MatPreallocated(B);
2203   PetscCheckSameComm(A,B);
2204   if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2205   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2206   if (A->M != B->M || A->N != B->N) SETERRQ4(PETSC_ERR_ARG_SIZ,"Mat A,Mat B: global dim %d %d",A->M,B->M,
2207                                              A->N,B->N);
2208 
2209   ierr = PetscLogEventBegin(MAT_Copy,A,B,0,0);CHKERRQ(ierr);
2210   if (A->ops->copy) {
2211     ierr = (*A->ops->copy)(A,B,str);CHKERRQ(ierr);
2212   } else { /* generic conversion */
2213     ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr);
2214   }
2215   ierr = PetscLogEventEnd(MAT_Copy,A,B,0,0);CHKERRQ(ierr);
2216   PetscFunctionReturn(0);
2217 }
2218 
2219 #include "petscsys.h"
2220 PetscTruth MatConvertRegisterAllCalled = PETSC_FALSE;
2221 PetscFList      MatConvertList              = 0;
2222 
2223 #undef __FUNC__
2224 #define __FUNC__ "MatConvertRegister"
2225 /*@C
2226     MatConvertRegister - Allows one to register a routine that reads matrices
2227         from a binary file for a particular matrix type.
2228 
2229   Not Collective
2230 
2231   Input Parameters:
2232 +   type - the type of matrix (defined in include/petscmat.h), for example, MATSEQAIJ.
2233 -   Converter - the function that reads the matrix from the binary file.
2234 
2235   Level: developer
2236 
2237 .seealso: MatConvertRegisterAll(), MatConvert()
2238 
2239 @*/
2240 int MatConvertRegister(char *sname,char *path,char *name,int (*function)(Mat,MatType,Mat*))
2241 {
2242   int  ierr;
2243   char fullname[256];
2244 
2245   PetscFunctionBegin;
2246   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
2247   ierr = PetscFListAdd(&MatConvertList,sname,fullname,(int (*)(void*))function);CHKERRQ(ierr);
2248   PetscFunctionReturn(0);
2249 }
2250 
2251 #undef __FUNC__
2252 #define __FUNC__ "MatConvert"
2253 /*@C
2254    MatConvert - Converts a matrix to another matrix, either of the same
2255    or different type.
2256 
2257    Collective on Mat
2258 
2259    Input Parameters:
2260 +  mat - the matrix
2261 -  newtype - new matrix type.  Use MATSAME to create a new matrix of the
2262    same type as the original matrix.
2263 
2264    Output Parameter:
2265 .  M - pointer to place new matrix
2266 
2267    Notes:
2268    MatConvert() first creates a new matrix and then copies the data from
2269    the first matrix.  A related routine is MatCopy(), which copies the matrix
2270    entries of one matrix to another already existing matrix context.
2271 
2272    Level: intermediate
2273 
2274    Concepts: matrices^converting between storage formats
2275 
2276 .seealso: MatCopy(), MatDuplicate()
2277 @*/
2278 int MatConvert(Mat mat,MatType newtype,Mat *M)
2279 {
2280   int        ierr;
2281   PetscTruth sametype,issame,flg;
2282   char       convname[256],mtype[256];
2283 
2284   PetscFunctionBegin;
2285   PetscValidHeaderSpecific(mat,MAT_COOKIE);
2286   PetscValidType(mat);
2287   MatPreallocated(mat);
2288   PetscValidPointer(M);
2289   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2290   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2291 
2292   ierr = PetscOptionsGetString(PETSC_NULL,"-matconvert_type",mtype,256,&flg);CHKERRQ(ierr);
2293   if (flg) {
2294     newtype = mtype;
2295   }
2296   ierr = PetscLogEventBegin(MAT_Convert,mat,0,0,0);CHKERRQ(ierr);
2297 
2298   ierr = PetscTypeCompare((PetscObject)mat,newtype,&sametype);CHKERRQ(ierr);
2299   ierr = PetscStrcmp(newtype,"same",&issame);CHKERRQ(ierr);
2300   if ((sametype || issame) && mat->ops->duplicate) {
2301     ierr = (*mat->ops->duplicate)(mat,MAT_COPY_VALUES,M);CHKERRQ(ierr);
2302   } else {
2303     int (*conv)(Mat,MatType,Mat*);
2304     ierr = PetscStrcpy(convname,"MatConvertTo_");CHKERRQ(ierr);
2305     ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr);
2306     ierr =  PetscFListFind(mat->comm,MatConvertList,convname,(int(**)(void*))&conv);CHKERRQ(ierr);
2307     if (conv) {
2308       ierr = (*conv)(mat,newtype,M);CHKERRQ(ierr);
2309     } else {
2310       ierr = PetscStrcpy(convname,"MatConvert_");CHKERRQ(ierr);
2311       ierr = PetscStrcat(convname,mat->type_name);CHKERRQ(ierr);
2312       ierr = PetscStrcat(convname,"_");CHKERRQ(ierr);
2313       ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr);
2314       ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr);
2315       ierr = PetscObjectQueryFunction((PetscObject)mat,convname,(void**)&conv);CHKERRQ(ierr);
2316       if (conv) {
2317         ierr = (*conv)(mat,newtype,M);CHKERRQ(ierr);
2318       } else {
2319         if (mat->ops->convert) {
2320           ierr = (*mat->ops->convert)(mat,newtype,M);CHKERRQ(ierr);
2321         } else {
2322           ierr = MatConvert_Basic(mat,newtype,M);CHKERRQ(ierr);
2323         }
2324       }
2325     }
2326   }
2327   ierr = PetscLogEventEnd(MAT_Convert,mat,0,0,0);CHKERRQ(ierr);
2328   PetscFunctionReturn(0);
2329 }
2330 
2331 
2332 #undef __FUNC__
2333 #define __FUNC__ "MatDuplicate"
2334 /*@C
2335    MatDuplicate - Duplicates a matrix including the non-zero structure.
2336 
2337    Collective on Mat
2338 
2339    Input Parameters:
2340 +  mat - the matrix
2341 -  op - either MAT_DO_NOT_COPY_VALUES or MAT_COPY_VALUES, cause it to copy nonzero
2342         values as well or not
2343 
2344    Output Parameter:
2345 .  M - pointer to place new matrix
2346 
2347    Level: intermediate
2348 
2349    Concepts: matrices^duplicating
2350 
2351 .seealso: MatCopy(), MatConvert()
2352 @*/
2353 int MatDuplicate(Mat mat,MatDuplicateOption op,Mat *M)
2354 {
2355   int ierr;
2356 
2357   PetscFunctionBegin;
2358   PetscValidHeaderSpecific(mat,MAT_COOKIE);
2359   PetscValidType(mat);
2360   MatPreallocated(mat);
2361   PetscValidPointer(M);
2362   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2363   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2364 
2365   *M  = 0;
2366   ierr = PetscLogEventBegin(MAT_Convert,mat,0,0,0);CHKERRQ(ierr);
2367   if (!mat->ops->duplicate) {
2368     SETERRQ(PETSC_ERR_SUP,"Not written for this matrix type");
2369   }
2370   ierr = (*mat->ops->duplicate)(mat,op,M);CHKERRQ(ierr);
2371   ierr = PetscLogEventEnd(MAT_Convert,mat,0,0,0);CHKERRQ(ierr);
2372   PetscFunctionReturn(0);
2373 }
2374 
2375 #undef __FUNC__
2376 #define __FUNC__ "MatGetDiagonal"
2377 /*@
2378    MatGetDiagonal - Gets the diagonal of a matrix.
2379 
2380    Collective on Mat and Vec
2381 
2382    Input Parameters:
2383 +  mat - the matrix
2384 -  v - the vector for storing the diagonal
2385 
2386    Output Parameter:
2387 .  v - the diagonal of the matrix
2388 
2389    Notes:
2390    For the SeqAIJ matrix format, this routine may also be called
2391    on a LU factored matrix; in that case it routines the reciprocal of
2392    the diagonal entries in U. It returns the entries permuted by the
2393    row and column permutation used during the symbolic factorization.
2394 
2395    Level: intermediate
2396 
2397    Concepts: matrices^accessing diagonals
2398 
2399 .seealso: MatGetRow(), MatGetSubmatrices(), MatGetSubmatrix(), MatGetRowMax()
2400 @*/
2401 int MatGetDiagonal(Mat mat,Vec v)
2402 {
2403   int ierr;
2404 
2405   PetscFunctionBegin;
2406   PetscValidHeaderSpecific(mat,MAT_COOKIE);
2407   PetscValidType(mat);
2408   MatPreallocated(mat);
2409   PetscValidHeaderSpecific(v,VEC_COOKIE);
2410   /* PetscCheckSameComm(mat,v); Could be MPI vector but Seq matrix cause of two submatrix storage */
2411   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2412   if (!mat->ops->getdiagonal) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name);
2413 
2414   ierr = (*mat->ops->getdiagonal)(mat,v);CHKERRQ(ierr);
2415   PetscFunctionReturn(0);
2416 }
2417 
2418 #undef __FUNC__
2419 #define __FUNC__ "MatGetRowMax"
2420 /*@
2421    MatGetRowMax - Gets the maximum value (in absolute value) of each
2422         row of the matrix
2423 
2424    Collective on Mat and Vec
2425 
2426    Input Parameters:
2427 .  mat - the matrix
2428 
2429    Output Parameter:
2430 .  v - the vector for storing the maximums
2431 
2432    Level: intermediate
2433 
2434    Concepts: matrices^getting row maximums
2435 
2436 .seealso: MatGetDiagonal(), MatGetSubmatrices(), MatGetSubmatrix()
2437 @*/
2438 int MatGetRowMax(Mat mat,Vec v)
2439 {
2440   int ierr;
2441 
2442   PetscFunctionBegin;
2443   PetscValidHeaderSpecific(mat,MAT_COOKIE);
2444   PetscValidType(mat);
2445   MatPreallocated(mat);
2446   PetscValidHeaderSpecific(v,VEC_COOKIE);
2447   /* PetscCheckSameComm(mat,v); Could be MPI vector but Seq matrix cause of two submatrix storage */
2448   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2449   if (!mat->ops->getrowmax) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name);
2450 
2451   ierr = (*mat->ops->getrowmax)(mat,v);CHKERRQ(ierr);
2452   PetscFunctionReturn(0);
2453 }
2454 
2455 #undef __FUNC__
2456 #define __FUNC__ "MatTranspose"
2457 /*@C
2458    MatTranspose - Computes an in-place or out-of-place transpose of a matrix.
2459 
2460    Collective on Mat
2461 
2462    Input Parameter:
2463 .  mat - the matrix to transpose
2464 
2465    Output Parameters:
2466 .  B - the transpose (or pass in PETSC_NULL for an in-place transpose)
2467 
2468    Level: intermediate
2469 
2470    Concepts: matrices^transposing
2471 
2472 .seealso: MatMultTranspose(), MatMultTransposeAdd()
2473 @*/
2474 int MatTranspose(Mat mat,Mat *B)
2475 {
2476   int ierr;
2477 
2478   PetscFunctionBegin;
2479   PetscValidHeaderSpecific(mat,MAT_COOKIE);
2480   PetscValidType(mat);
2481   MatPreallocated(mat);
2482   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2483   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2484   if (!mat->ops->transpose) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name);
2485   ierr = (*mat->ops->transpose)(mat,B);CHKERRQ(ierr);
2486   PetscFunctionReturn(0);
2487 }
2488 
2489 #undef __FUNC__
2490 #define __FUNC__ "MatPermute"
2491 /*@C
2492    MatPermute - Creates a new matrix with rows and columns permuted from the
2493    original.
2494 
2495    Collective on Mat
2496 
2497    Input Parameters:
2498 +  mat - the matrix to permute
2499 .  row - row permutation, each processor supplies only the permutation for its rows
2500 -  col - column permutation, each processor needs the entire column permutation, that is
2501          this is the same size as the total number of columns in the matrix
2502 
2503    Output Parameters:
2504 .  B - the permuted matrix
2505 
2506    Level: advanced
2507 
2508    Concepts: matrices^permuting
2509 
2510 .seealso: MatGetOrdering()
2511 @*/
2512 int MatPermute(Mat mat,IS row,IS col,Mat *B)
2513 {
2514   int ierr;
2515 
2516   PetscFunctionBegin;
2517   PetscValidHeaderSpecific(mat,MAT_COOKIE);
2518   PetscValidType(mat);
2519   MatPreallocated(mat);
2520   PetscValidHeaderSpecific(row,IS_COOKIE);
2521   PetscValidHeaderSpecific(col,IS_COOKIE);
2522   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2523   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2524   if (!mat->ops->permute) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name);
2525   ierr = (*mat->ops->permute)(mat,row,col,B);CHKERRQ(ierr);
2526   PetscFunctionReturn(0);
2527 }
2528 
2529 #undef __FUNC__
2530 #define __FUNC__ "MatEqual"
2531 /*@
2532    MatEqual - Compares two matrices.
2533 
2534    Collective on Mat
2535 
2536    Input Parameters:
2537 +  A - the first matrix
2538 -  B - the second matrix
2539 
2540    Output Parameter:
2541 .  flg - PETSC_TRUE if the matrices are equal; PETSC_FALSE otherwise.
2542 
2543    Level: intermediate
2544 
2545    Concepts: matrices^equality between
2546 @*/
2547 int MatEqual(Mat A,Mat B,PetscTruth *flg)
2548 {
2549   int ierr;
2550 
2551   PetscFunctionBegin;
2552   PetscValidHeaderSpecific(A,MAT_COOKIE);
2553   PetscValidHeaderSpecific(B,MAT_COOKIE);
2554   PetscValidType(A);
2555   MatPreallocated(A);
2556   PetscValidType(B);
2557   MatPreallocated(B);
2558   PetscValidIntPointer(flg);
2559   PetscCheckSameComm(A,B);
2560   if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2561   if (!B->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2562   if (A->M != B->M || A->N != B->N) SETERRQ4(PETSC_ERR_ARG_SIZ,"Mat A,Mat B: global dim %d %d %d %d",A->M,B->M,A->N,B->N);
2563   if (!A->ops->equal) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",A->type_name);
2564   ierr = (*A->ops->equal)(A,B,flg);CHKERRQ(ierr);
2565   PetscFunctionReturn(0);
2566 }
2567 
2568 #undef __FUNC__
2569 #define __FUNC__ "MatDiagonalScale"
2570 /*@
2571    MatDiagonalScale - Scales a matrix on the left and right by diagonal
2572    matrices that are stored as vectors.  Either of the two scaling
2573    matrices can be PETSC_NULL.
2574 
2575    Collective on Mat
2576 
2577    Input Parameters:
2578 +  mat - the matrix to be scaled
2579 .  l - the left scaling vector (or PETSC_NULL)
2580 -  r - the right scaling vector (or PETSC_NULL)
2581 
2582    Notes:
2583    MatDiagonalScale() computes A = LAR, where
2584    L = a diagonal matrix, R = a diagonal matrix
2585 
2586    Level: intermediate
2587 
2588    Concepts: matrices^diagonal scaling
2589    Concepts: diagonal scaling of matrices
2590 
2591 .seealso: MatScale()
2592 @*/
2593 int MatDiagonalScale(Mat mat,Vec l,Vec r)
2594 {
2595   int ierr;
2596 
2597   PetscFunctionBegin;
2598   PetscValidHeaderSpecific(mat,MAT_COOKIE);
2599   PetscValidType(mat);
2600   MatPreallocated(mat);
2601   if (!mat->ops->diagonalscale) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name);
2602   if (l) {PetscValidHeaderSpecific(l,VEC_COOKIE);PetscCheckSameComm(mat,l);}
2603   if (r) {PetscValidHeaderSpecific(r,VEC_COOKIE);PetscCheckSameComm(mat,r);}
2604   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2605   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2606 
2607   ierr = PetscLogEventBegin(MAT_Scale,mat,0,0,0);CHKERRQ(ierr);
2608   ierr = (*mat->ops->diagonalscale)(mat,l,r);CHKERRQ(ierr);
2609   ierr = PetscLogEventEnd(MAT_Scale,mat,0,0,0);CHKERRQ(ierr);
2610   PetscFunctionReturn(0);
2611 }
2612 
2613 #undef __FUNC__
2614 #define __FUNC__ "MatScale"
2615 /*@
2616     MatScale - Scales all elements of a matrix by a given number.
2617 
2618     Collective on Mat
2619 
2620     Input Parameters:
2621 +   mat - the matrix to be scaled
2622 -   a  - the scaling value
2623 
2624     Output Parameter:
2625 .   mat - the scaled matrix
2626 
2627     Level: intermediate
2628 
2629     Concepts: matrices^scaling all entries
2630 
2631 .seealso: MatDiagonalScale()
2632 @*/
2633 int MatScale(Scalar *a,Mat mat)
2634 {
2635   int ierr;
2636 
2637   PetscFunctionBegin;
2638   PetscValidHeaderSpecific(mat,MAT_COOKIE);
2639   PetscValidType(mat);
2640   MatPreallocated(mat);
2641   PetscValidScalarPointer(a);
2642   if (!mat->ops->scale) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name);
2643   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2644   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2645 
2646   ierr = PetscLogEventBegin(MAT_Scale,mat,0,0,0);CHKERRQ(ierr);
2647   ierr = (*mat->ops->scale)(a,mat);CHKERRQ(ierr);
2648   ierr = PetscLogEventEnd(MAT_Scale,mat,0,0,0);CHKERRQ(ierr);
2649   PetscFunctionReturn(0);
2650 }
2651 
2652 #undef __FUNC__
2653 #define __FUNC__ "MatNorm"
2654 /*@
2655    MatNorm - Calculates various norms of a matrix.
2656 
2657    Collective on Mat
2658 
2659    Input Parameters:
2660 +  mat - the matrix
2661 -  type - the type of norm, NORM_1, NORM_2, NORM_FROBENIUS, NORM_INFINITY
2662 
2663    Output Parameters:
2664 .  norm - the resulting norm
2665 
2666    Level: intermediate
2667 
2668    Concepts: matrices^norm
2669    Concepts: norm^of matrix
2670 @*/
2671 int MatNorm(Mat mat,NormType type,PetscReal *norm)
2672 {
2673   int ierr;
2674 
2675   PetscFunctionBegin;
2676   PetscValidHeaderSpecific(mat,MAT_COOKIE);
2677   PetscValidType(mat);
2678   MatPreallocated(mat);
2679   PetscValidScalarPointer(norm);
2680 
2681   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2682   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2683   if (!mat->ops->norm) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name);
2684   ierr = (*mat->ops->norm)(mat,type,norm);CHKERRQ(ierr);
2685   PetscFunctionReturn(0);
2686 }
2687 
2688 /*
2689      This variable is used to prevent counting of MatAssemblyBegin() that
2690    are called from within a MatAssemblyEnd().
2691 */
2692 static int MatAssemblyEnd_InUse = 0;
2693 #undef __FUNC__
2694 #define __FUNC__ "MatAssemblyBegin"
2695 /*@
2696    MatAssemblyBegin - Begins assembling the matrix.  This routine should
2697    be called after completing all calls to MatSetValues().
2698 
2699    Collective on Mat
2700 
2701    Input Parameters:
2702 +  mat - the matrix
2703 -  type - type of assembly, either MAT_FLUSH_ASSEMBLY or MAT_FINAL_ASSEMBLY
2704 
2705    Notes:
2706    MatSetValues() generally caches the values.  The matrix is ready to
2707    use only after MatAssemblyBegin() and MatAssemblyEnd() have been called.
2708    Use MAT_FLUSH_ASSEMBLY when switching between ADD_VALUES and INSERT_VALUES
2709    in MatSetValues(); use MAT_FINAL_ASSEMBLY for the final assembly before
2710    using the matrix.
2711 
2712    Level: beginner
2713 
2714    Concepts: matrices^assembling
2715 
2716 .seealso: MatAssemblyEnd(), MatSetValues(), MatAssembled()
2717 @*/
2718 int MatAssemblyBegin(Mat mat,MatAssemblyType type)
2719 {
2720   int ierr;
2721 
2722   PetscFunctionBegin;
2723   PetscValidHeaderSpecific(mat,MAT_COOKIE);
2724   PetscValidType(mat);
2725   MatPreallocated(mat);
2726   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix.\nDid you forget to call MatSetUnfactored()?");
2727   if (mat->assembled) {
2728     mat->was_assembled = PETSC_TRUE;
2729     mat->assembled     = PETSC_FALSE;
2730   }
2731   if (!MatAssemblyEnd_InUse) {
2732     ierr = PetscLogEventBegin(MAT_AssemblyBegin,mat,0,0,0);CHKERRQ(ierr);
2733     if (mat->ops->assemblybegin){ierr = (*mat->ops->assemblybegin)(mat,type);CHKERRQ(ierr);}
2734     ierr = PetscLogEventEnd(MAT_AssemblyBegin,mat,0,0,0);CHKERRQ(ierr);
2735   } else {
2736     if (mat->ops->assemblybegin){ierr = (*mat->ops->assemblybegin)(mat,type);CHKERRQ(ierr);}
2737   }
2738   PetscFunctionReturn(0);
2739 }
2740 
2741 #undef __FUNC__
2742 #define __FUNC__ "MatAssembed"
2743 /*@
2744    MatAssembled - Indicates if a matrix has been assembled and is ready for
2745      use; for example, in matrix-vector product.
2746 
2747    Collective on Mat
2748 
2749    Input Parameter:
2750 .  mat - the matrix
2751 
2752    Output Parameter:
2753 .  assembled - PETSC_TRUE or PETSC_FALSE
2754 
2755    Level: advanced
2756 
2757    Concepts: matrices^assembled?
2758 
2759 .seealso: MatAssemblyEnd(), MatSetValues(), MatAssemblyBegin()
2760 @*/
2761 int MatAssembled(Mat mat,PetscTruth *assembled)
2762 {
2763   PetscFunctionBegin;
2764   PetscValidHeaderSpecific(mat,MAT_COOKIE);
2765   PetscValidType(mat);
2766   MatPreallocated(mat);
2767   *assembled = mat->assembled;
2768   PetscFunctionReturn(0);
2769 }
2770 
2771 #undef __FUNC__
2772 #define __FUNC__ "MatView_Private"
2773 /*
2774     Processes command line options to determine if/how a matrix
2775   is to be viewed. Called by MatAssemblyEnd() and MatLoad().
2776 */
2777 int MatView_Private(Mat mat)
2778 {
2779   int        ierr;
2780   PetscTruth flg;
2781 
2782   PetscFunctionBegin;
2783   ierr = PetscOptionsHasName(mat->prefix,"-mat_view_info",&flg);CHKERRQ(ierr);
2784   if (flg) {
2785     ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_(mat->comm),PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
2786     ierr = MatView(mat,PETSC_VIEWER_STDOUT_(mat->comm));CHKERRQ(ierr);
2787     ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_(mat->comm));CHKERRQ(ierr);
2788   }
2789   ierr = PetscOptionsHasName(mat->prefix,"-mat_view_info_detailed",&flg);CHKERRQ(ierr);
2790   if (flg) {
2791     ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_(mat->comm),PETSC_VIEWER_ASCII_INFO_LONG);CHKERRQ(ierr);
2792     ierr = MatView(mat,PETSC_VIEWER_STDOUT_(mat->comm));CHKERRQ(ierr);
2793     ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_(mat->comm));CHKERRQ(ierr);
2794   }
2795   ierr = PetscOptionsHasName(mat->prefix,"-mat_view",&flg);CHKERRQ(ierr);
2796   if (flg) {
2797     ierr = MatView(mat,PETSC_VIEWER_STDOUT_(mat->comm));CHKERRQ(ierr);
2798   }
2799   ierr = PetscOptionsHasName(mat->prefix,"-mat_view_matlab",&flg);CHKERRQ(ierr);
2800   if (flg) {
2801     ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_(mat->comm),PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
2802     ierr = MatView(mat,PETSC_VIEWER_STDOUT_(mat->comm));CHKERRQ(ierr);
2803     ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_(mat->comm));CHKERRQ(ierr);
2804   }
2805   ierr = PetscOptionsHasName(mat->prefix,"-mat_view_draw",&flg);CHKERRQ(ierr);
2806   if (flg) {
2807     ierr = PetscOptionsHasName(mat->prefix,"-mat_view_contour",&flg);CHKERRQ(ierr);
2808     if (flg) {
2809       PetscViewerPushFormat(PETSC_VIEWER_DRAW_(mat->comm),PETSC_VIEWER_DRAW_CONTOUR);CHKERRQ(ierr);
2810     }
2811     ierr = MatView(mat,PETSC_VIEWER_DRAW_(mat->comm));CHKERRQ(ierr);
2812     ierr = PetscViewerFlush(PETSC_VIEWER_DRAW_(mat->comm));CHKERRQ(ierr);
2813     if (flg) {
2814       PetscViewerPopFormat(PETSC_VIEWER_DRAW_(mat->comm));CHKERRQ(ierr);
2815     }
2816   }
2817   ierr = PetscOptionsHasName(mat->prefix,"-mat_view_socket",&flg);CHKERRQ(ierr);
2818   if (flg) {
2819     ierr = MatView(mat,PETSC_VIEWER_SOCKET_(mat->comm));CHKERRQ(ierr);
2820     ierr = PetscViewerFlush(PETSC_VIEWER_SOCKET_(mat->comm));CHKERRQ(ierr);
2821   }
2822   ierr = PetscOptionsHasName(mat->prefix,"-mat_view_binary",&flg);CHKERRQ(ierr);
2823   if (flg) {
2824     ierr = MatView(mat,PETSC_VIEWER_BINARY_(mat->comm));CHKERRQ(ierr);
2825     ierr = PetscViewerFlush(PETSC_VIEWER_BINARY_(mat->comm));CHKERRQ(ierr);
2826   }
2827   PetscFunctionReturn(0);
2828 }
2829 
2830 #undef __FUNC__
2831 #define __FUNC__ "MatAssemblyEnd"
2832 /*@
2833    MatAssemblyEnd - Completes assembling the matrix.  This routine should
2834    be called after MatAssemblyBegin().
2835 
2836    Collective on Mat
2837 
2838    Input Parameters:
2839 +  mat - the matrix
2840 -  type - type of assembly, either MAT_FLUSH_ASSEMBLY or MAT_FINAL_ASSEMBLY
2841 
2842    Options Database Keys:
2843 +  -mat_view_info - Prints info on matrix at conclusion of MatEndAssembly()
2844 .  -mat_view_info_detailed - Prints more detailed info
2845 .  -mat_view - Prints matrix in ASCII format
2846 .  -mat_view_matlab - Prints matrix in Matlab format
2847 .  -mat_view_draw - PetscDraws nonzero structure of matrix, using MatView() and PetscDrawOpenX().
2848 .  -display <name> - Sets display name (default is host)
2849 -  -draw_pause <sec> - Sets number of seconds to pause after display
2850 
2851    Notes:
2852    MatSetValues() generally caches the values.  The matrix is ready to
2853    use only after MatAssemblyBegin() and MatAssemblyEnd() have been called.
2854    Use MAT_FLUSH_ASSEMBLY when switching between ADD_VALUES and INSERT_VALUES
2855    in MatSetValues(); use MAT_FINAL_ASSEMBLY for the final assembly before
2856    using the matrix.
2857 
2858    Level: beginner
2859 
2860 .seealso: MatAssemblyBegin(), MatSetValues(), PetscDrawOpenX(), MatView(), MatAssembled()
2861 @*/
2862 int MatAssemblyEnd(Mat mat,MatAssemblyType type)
2863 {
2864   int        ierr;
2865   static int inassm = 0;
2866 
2867   PetscFunctionBegin;
2868   PetscValidHeaderSpecific(mat,MAT_COOKIE);
2869   PetscValidType(mat);
2870   MatPreallocated(mat);
2871 
2872   inassm++;
2873   MatAssemblyEnd_InUse++;
2874   if (MatAssemblyEnd_InUse == 1) { /* Do the logging only the first time through */
2875     ierr = PetscLogEventBegin(MAT_AssemblyEnd,mat,0,0,0);CHKERRQ(ierr);
2876     if (mat->ops->assemblyend) {
2877       ierr = (*mat->ops->assemblyend)(mat,type);CHKERRQ(ierr);
2878     }
2879     ierr = PetscLogEventEnd(MAT_AssemblyEnd,mat,0,0,0);CHKERRQ(ierr);
2880   } else {
2881     if (mat->ops->assemblyend) {
2882       ierr = (*mat->ops->assemblyend)(mat,type);CHKERRQ(ierr);
2883     }
2884   }
2885 
2886   /* Flush assembly is not a true assembly */
2887   if (type != MAT_FLUSH_ASSEMBLY) {
2888     mat->assembled  = PETSC_TRUE; mat->num_ass++;
2889   }
2890   mat->insertmode = NOT_SET_VALUES;
2891   MatAssemblyEnd_InUse--;
2892 
2893   if (inassm == 1 && type != MAT_FLUSH_ASSEMBLY) {
2894     ierr = MatView_Private(mat);CHKERRQ(ierr);
2895   }
2896   inassm--;
2897   PetscFunctionReturn(0);
2898 }
2899 
2900 
2901 #undef __FUNC__
2902 #define __FUNC__ "MatCompress"
2903 /*@
2904    MatCompress - Tries to store the matrix in as little space as
2905    possible.  May fail if memory is already fully used, since it
2906    tries to allocate new space.
2907 
2908    Collective on Mat
2909 
2910    Input Parameters:
2911 .  mat - the matrix
2912 
2913    Level: advanced
2914 
2915 @*/
2916 int MatCompress(Mat mat)
2917 {
2918   int ierr;
2919 
2920   PetscFunctionBegin;
2921   PetscValidHeaderSpecific(mat,MAT_COOKIE);
2922   PetscValidType(mat);
2923   MatPreallocated(mat);
2924   if (mat->ops->compress) {ierr = (*mat->ops->compress)(mat);CHKERRQ(ierr);}
2925   PetscFunctionReturn(0);
2926 }
2927 
2928 #undef __FUNC__
2929 #define __FUNC__ "MatSetOption"
2930 /*@
2931    MatSetOption - Sets a parameter option for a matrix. Some options
2932    may be specific to certain storage formats.  Some options
2933    determine how values will be inserted (or added). Sorted,
2934    row-oriented input will generally assemble the fastest. The default
2935    is row-oriented, nonsorted input.
2936 
2937    Collective on Mat
2938 
2939    Input Parameters:
2940 +  mat - the matrix
2941 -  option - the option, one of those listed below (and possibly others),
2942              e.g., MAT_ROWS_SORTED, MAT_NEW_NONZERO_LOCATION_ERR
2943 
2944    Options Describing Matrix Structure:
2945 +    MAT_SYMMETRIC - symmetric in terms of both structure and value
2946 -    MAT_STRUCTURALLY_SYMMETRIC - symmetric nonzero structure
2947 
2948    Options For Use with MatSetValues():
2949    Insert a logically dense subblock, which can be
2950 +    MAT_ROW_ORIENTED - row-oriented
2951 .    MAT_COLUMN_ORIENTED - column-oriented
2952 .    MAT_ROWS_SORTED - sorted by row
2953 .    MAT_ROWS_UNSORTED - not sorted by row
2954 .    MAT_COLUMNS_SORTED - sorted by column
2955 -    MAT_COLUMNS_UNSORTED - not sorted by column
2956 
2957    Not these options reflect the data you pass in with MatSetValues(); it has
2958    nothing to do with how the data is stored internally in the matrix
2959    data structure.
2960 
2961    When (re)assembling a matrix, we can restrict the input for
2962    efficiency/debugging purposes.  These options include
2963 +    MAT_NO_NEW_NONZERO_LOCATIONS - additional insertions will not be
2964         allowed if they generate a new nonzero
2965 .    MAT_YES_NEW_NONZERO_LOCATIONS - additional insertions will be allowed
2966 .    MAT_NO_NEW_DIAGONALS - additional insertions will not be allowed if
2967          they generate a nonzero in a new diagonal (for block diagonal format only)
2968 .    MAT_YES_NEW_DIAGONALS - new diagonals will be allowed (for block diagonal format only)
2969 .    MAT_IGNORE_OFF_PROC_ENTRIES - drops off-processor entries
2970 .    MAT_NEW_NONZERO_LOCATION_ERR - generates an error for new matrix entry
2971 -    MAT_USE_HASH_TABLE - uses a hash table to speed up matrix assembly
2972 
2973    Notes:
2974    Some options are relevant only for particular matrix types and
2975    are thus ignored by others.  Other options are not supported by
2976    certain matrix types and will generate an error message if set.
2977 
2978    If using a Fortran 77 module to compute a matrix, one may need to
2979    use the column-oriented option (or convert to the row-oriented
2980    format).
2981 
2982    MAT_NO_NEW_NONZERO_LOCATIONS indicates that any add or insertion
2983    that would generate a new entry in the nonzero structure is instead
2984    ignored.  Thus, if memory has not alredy been allocated for this particular
2985    data, then the insertion is ignored. For dense matrices, in which
2986    the entire array is allocated, no entries are ever ignored.
2987 
2988    MAT_NEW_NONZERO_LOCATION_ERR indicates that any add or insertion
2989    that would generate a new entry in the nonzero structure instead produces
2990    an error. (Currently supported for AIJ and BAIJ formats only.)
2991    This is a useful flag when using SAME_NONZERO_PATTERN in calling
2992    SLESSetOperators() to ensure that the nonzero pattern truely does
2993    remain unchanged.
2994 
2995    MAT_NEW_NONZERO_ALLOCATION_ERR indicates that any add or insertion
2996    that would generate a new entry that has not been preallocated will
2997    instead produce an error. (Currently supported for AIJ and BAIJ formats
2998    only.) This is a useful flag when debugging matrix memory preallocation.
2999 
3000    MAT_IGNORE_OFF_PROC_ENTRIES indicates entries destined for
3001    other processors should be dropped, rather than stashed.
3002    This is useful if you know that the "owning" processor is also
3003    always generating the correct matrix entries, so that PETSc need
3004    not transfer duplicate entries generated on another processor.
3005 
3006    MAT_USE_HASH_TABLE indicates that a hash table be used to improve the
3007    searches during matrix assembly. When this flag is set, the hash table
3008    is created during the first Matrix Assembly. This hash table is
3009    used the next time through, during MatSetVaules()/MatSetVaulesBlocked()
3010    to improve the searching of indices. MAT_NO_NEW_NONZERO_LOCATIONS flag
3011    should be used with MAT_USE_HASH_TABLE flag. This option is currently
3012    supported by MATMPIBAIJ format only.
3013 
3014    MAT_KEEP_ZEROED_ROWS indicates when MatZeroRows() is called the zeroed entries
3015    are kept in the nonzero structure
3016 
3017    MAT_IGNORE_ZERO_ENTRIES - when using ADD_VALUES for AIJ matrices this will stop
3018    zero values from creating a zero location in the matrix
3019 
3020    MAT_USE_INODES - indicates using inode version of the code - works with AIJ and
3021    ROWBS matrix types
3022 
3023    MAT_DO_NOT_USE_INODES - indicates not using inode version of the code - works
3024    with AIJ and ROWBS matrix types
3025 
3026    Level: intermediate
3027 
3028    Concepts: matrices^setting options
3029 
3030 @*/
3031 int MatSetOption(Mat mat,MatOption op)
3032 {
3033   int ierr;
3034 
3035   PetscFunctionBegin;
3036   PetscValidHeaderSpecific(mat,MAT_COOKIE);
3037   PetscValidType(mat);
3038   MatPreallocated(mat);
3039   if (mat->ops->setoption) {ierr = (*mat->ops->setoption)(mat,op);CHKERRQ(ierr);}
3040   PetscFunctionReturn(0);
3041 }
3042 
3043 #undef __FUNC__
3044 #define __FUNC__ "MatZeroEntries"
3045 /*@
3046    MatZeroEntries - Zeros all entries of a matrix.  For sparse matrices
3047    this routine retains the old nonzero structure.
3048 
3049    Collective on Mat
3050 
3051    Input Parameters:
3052 .  mat - the matrix
3053 
3054    Level: intermediate
3055 
3056    Concepts: matrices^zeroing
3057 
3058 .seealso: MatZeroRows()
3059 @*/
3060 int MatZeroEntries(Mat mat)
3061 {
3062   int ierr;
3063 
3064   PetscFunctionBegin;
3065   PetscValidHeaderSpecific(mat,MAT_COOKIE);
3066   PetscValidType(mat);
3067   MatPreallocated(mat);
3068   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
3069   if (!mat->ops->zeroentries) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name);
3070 
3071   ierr = PetscLogEventBegin(MAT_ZeroEntries,mat,0,0,0);CHKERRQ(ierr);
3072   ierr = (*mat->ops->zeroentries)(mat);CHKERRQ(ierr);
3073   ierr = PetscLogEventEnd(MAT_ZeroEntries,mat,0,0,0);CHKERRQ(ierr);
3074   PetscFunctionReturn(0);
3075 }
3076 
3077 #undef __FUNC__
3078 #define __FUNC__ "MatZeroRows"
3079 /*@C
3080    MatZeroRows - Zeros all entries (except possibly the main diagonal)
3081    of a set of rows of a matrix.
3082 
3083    Collective on Mat
3084 
3085    Input Parameters:
3086 +  mat - the matrix
3087 .  is - index set of rows to remove
3088 -  diag - pointer to value put in all diagonals of eliminated rows.
3089           Note that diag is not a pointer to an array, but merely a
3090           pointer to a single value.
3091 
3092    Notes:
3093    For the AIJ and BAIJ matrix formats this removes the old nonzero structure,
3094    but does not release memory.  For the dense and block diagonal
3095    formats this does not alter the nonzero structure.
3096 
3097    If the option MatSetOption(mat,MAT_KEEP_ZEROED_ROWS) the nonzero structure
3098    of the matrix is not changed (even for AIJ and BAIJ matrices) the values are
3099    merely zeroed.
3100 
3101    The user can set a value in the diagonal entry (or for the AIJ and
3102    row formats can optionally remove the main diagonal entry from the
3103    nonzero structure as well, by passing a null pointer (PETSC_NULL
3104    in C or PETSC_NULL_SCALAR in Fortran) as the final argument).
3105 
3106    For the parallel case, all processes that share the matrix (i.e.,
3107    those in the communicator used for matrix creation) MUST call this
3108    routine, regardless of whether any rows being zeroed are owned by
3109    them.
3110 
3111 
3112    Level: intermediate
3113 
3114    Concepts: matrices^zeroing rows
3115 
3116 .seealso: MatZeroEntries(), MatZeroRowsLocal(), MatSetOption()
3117 @*/
3118 int MatZeroRows(Mat mat,IS is,Scalar *diag)
3119 {
3120   int ierr;
3121 
3122   PetscFunctionBegin;
3123   PetscValidHeaderSpecific(mat,MAT_COOKIE);
3124   PetscValidType(mat);
3125   MatPreallocated(mat);
3126   PetscValidHeaderSpecific(is,IS_COOKIE);
3127   if (diag) PetscValidScalarPointer(diag);
3128   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
3129   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
3130   if (!mat->ops->zerorows) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name);
3131 
3132   ierr = (*mat->ops->zerorows)(mat,is,diag);CHKERRQ(ierr);
3133   ierr = MatView_Private(mat);CHKERRQ(ierr);
3134   PetscFunctionReturn(0);
3135 }
3136 
3137 #undef __FUNC__
3138 #define __FUNC__ "MatZeroRowsLocal"
3139 /*@C
3140    MatZeroRowsLocal - Zeros all entries (except possibly the main diagonal)
3141    of a set of rows of a matrix; using local numbering of rows.
3142 
3143    Collective on Mat
3144 
3145    Input Parameters:
3146 +  mat - the matrix
3147 .  is - index set of rows to remove
3148 -  diag - pointer to value put in all diagonals of eliminated rows.
3149           Note that diag is not a pointer to an array, but merely a
3150           pointer to a single value.
3151 
3152    Notes:
3153    Before calling MatZeroRowsLocal(), the user must first set the
3154    local-to-global mapping by calling MatSetLocalToGlobalMapping().
3155 
3156    For the AIJ matrix formats this removes the old nonzero structure,
3157    but does not release memory.  For the dense and block diagonal
3158    formats this does not alter the nonzero structure.
3159 
3160    The user can set a value in the diagonal entry (or for the AIJ and
3161    row formats can optionally remove the main diagonal entry from the
3162    nonzero structure as well, by passing a null pointer (PETSC_NULL
3163    in C or PETSC_NULL_SCALAR in Fortran) as the final argument).
3164 
3165    Level: intermediate
3166 
3167    Concepts: matrices^zeroing
3168 
3169 .seealso: MatZeroEntries(), MatZeroRows(), MatSetLocalToGlobalMapping
3170 @*/
3171 int MatZeroRowsLocal(Mat mat,IS is,Scalar *diag)
3172 {
3173   int ierr;
3174   IS  newis;
3175 
3176   PetscFunctionBegin;
3177   PetscValidHeaderSpecific(mat,MAT_COOKIE);
3178   PetscValidType(mat);
3179   MatPreallocated(mat);
3180   PetscValidHeaderSpecific(is,IS_COOKIE);
3181   if (diag) PetscValidScalarPointer(diag);
3182   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
3183   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
3184 
3185   if (mat->ops->zerorowslocal) {
3186     ierr = (*mat->ops->zerorowslocal)(mat,is,diag);CHKERRQ(ierr);
3187   } else {
3188     if (!mat->mapping) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Need to provide local to global mapping to matrix first");
3189     ierr = ISLocalToGlobalMappingApplyIS(mat->mapping,is,&newis);CHKERRQ(ierr);
3190     ierr = (*mat->ops->zerorows)(mat,newis,diag);CHKERRQ(ierr);
3191     ierr = ISDestroy(newis);CHKERRQ(ierr);
3192   }
3193   PetscFunctionReturn(0);
3194 }
3195 
3196 #undef __FUNC__
3197 #define __FUNC__ "MatGetSize"
3198 /*@
3199    MatGetSize - Returns the numbers of rows and columns in a matrix.
3200 
3201    Not Collective
3202 
3203    Input Parameter:
3204 .  mat - the matrix
3205 
3206    Output Parameters:
3207 +  m - the number of global rows
3208 -  n - the number of global columns
3209 
3210    Level: beginner
3211 
3212    Concepts: matrices^size
3213 
3214 .seealso: MatGetLocalSize()
3215 @*/
3216 int MatGetSize(Mat mat,int *m,int* n)
3217 {
3218   PetscFunctionBegin;
3219   PetscValidHeaderSpecific(mat,MAT_COOKIE);
3220   if (m) *m = mat->M;
3221   if (n) *n = mat->N;
3222   PetscFunctionReturn(0);
3223 }
3224 
3225 #undef __FUNC__
3226 #define __FUNC__ "MatGetLocalSize"
3227 /*@
3228    MatGetLocalSize - Returns the number of rows and columns in a matrix
3229    stored locally.  This information may be implementation dependent, so
3230    use with care.
3231 
3232    Not Collective
3233 
3234    Input Parameters:
3235 .  mat - the matrix
3236 
3237    Output Parameters:
3238 +  m - the number of local rows
3239 -  n - the number of local columns
3240 
3241    Level: beginner
3242 
3243    Concepts: matrices^local size
3244 
3245 .seealso: MatGetSize()
3246 @*/
3247 int MatGetLocalSize(Mat mat,int *m,int* n)
3248 {
3249   PetscFunctionBegin;
3250   PetscValidHeaderSpecific(mat,MAT_COOKIE);
3251   if (m) *m = mat->m;
3252   if (n) *n = mat->n;
3253   PetscFunctionReturn(0);
3254 }
3255 
3256 #undef __FUNC__
3257 #define __FUNC__ "MatGetOwnershipRange"
3258 /*@
3259    MatGetOwnershipRange - Returns the range of matrix rows owned by
3260    this processor, assuming that the matrix is laid out with the first
3261    n1 rows on the first processor, the next n2 rows on the second, etc.
3262    For certain parallel layouts this range may not be well defined.
3263 
3264    Not Collective
3265 
3266    Input Parameters:
3267 .  mat - the matrix
3268 
3269    Output Parameters:
3270 +  m - the global index of the first local row
3271 -  n - one more than the global index of the last local row
3272 
3273    Level: beginner
3274 
3275    Concepts: matrices^row ownership
3276 @*/
3277 int MatGetOwnershipRange(Mat mat,int *m,int* n)
3278 {
3279   int ierr;
3280 
3281   PetscFunctionBegin;
3282   PetscValidHeaderSpecific(mat,MAT_COOKIE);
3283   PetscValidType(mat);
3284   MatPreallocated(mat);
3285   if (m) PetscValidIntPointer(m);
3286   if (n) PetscValidIntPointer(n);
3287   if (!mat->ops->getownershiprange) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name);
3288   ierr = (*mat->ops->getownershiprange)(mat,m,n);CHKERRQ(ierr);
3289   PetscFunctionReturn(0);
3290 }
3291 
3292 #undef __FUNC__
3293 #define __FUNC__ "MatILUFactorSymbolic"
3294 /*@
3295    MatILUFactorSymbolic - Performs symbolic ILU factorization of a matrix.
3296    Uses levels of fill only, not drop tolerance. Use MatLUFactorNumeric()
3297    to complete the factorization.
3298 
3299    Collective on Mat
3300 
3301    Input Parameters:
3302 +  mat - the matrix
3303 .  row - row permutation
3304 .  column - column permutation
3305 -  info - structure containing
3306 $      levels - number of levels of fill.
3307 $      expected fill - as ratio of original fill.
3308 $      1 or 0 - indicating force fill on diagonal (improves robustness for matrices
3309                 missing diagonal entries)
3310 
3311    Output Parameters:
3312 .  fact - new matrix that has been symbolically factored
3313 
3314    Notes:
3315    See the users manual for additional information about
3316    choosing the fill factor for better efficiency.
3317 
3318    Most users should employ the simplified SLES interface for linear solvers
3319    instead of working directly with matrix algebra routines such as this.
3320    See, e.g., SLESCreate().
3321 
3322    Level: developer
3323 
3324   Concepts: matrices^symbolic LU factorization
3325   Concepts: matrices^factorization
3326   Concepts: LU^symbolic factorization
3327 
3328 .seealso: MatLUFactorSymbolic(), MatLUFactorNumeric(), MatCholeskyFactor()
3329           MatGetOrdering()
3330 
3331 @*/
3332 int MatILUFactorSymbolic(Mat mat,IS row,IS col,MatILUInfo *info,Mat *fact)
3333 {
3334   int ierr;
3335 
3336   PetscFunctionBegin;
3337   PetscValidHeaderSpecific(mat,MAT_COOKIE);
3338   PetscValidType(mat);
3339   MatPreallocated(mat);
3340   PetscValidPointer(fact);
3341   if (info && info->levels < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Levels of fill negative %d",info->levels);
3342   if (info && info->fill < 1.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Expected fill less than 1.0 %g",info->fill);
3343   if (!mat->ops->ilufactorsymbolic) SETERRQ1(PETSC_ERR_SUP,"Matrix type %s  symbolic ILU",mat->type_name);
3344   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
3345   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
3346 
3347   ierr = PetscLogEventBegin(MAT_ILUFactorSymbolic,mat,row,col,0);CHKERRQ(ierr);
3348   ierr = (*mat->ops->ilufactorsymbolic)(mat,row,col,info,fact);CHKERRQ(ierr);
3349   ierr = PetscLogEventEnd(MAT_ILUFactorSymbolic,mat,row,col,0);CHKERRQ(ierr);
3350   PetscFunctionReturn(0);
3351 }
3352 
3353 #undef __FUNC__
3354 #define __FUNC__ "MatIncompleteCholeskyFactorSymbolic"
3355 /*@
3356    MatIncompleteCholeskyFactorSymbolic - Performs symbolic incomplete
3357    Cholesky factorization for a symmetric matrix.  Use
3358    MatCholeskyFactorNumeric() to complete the factorization.
3359 
3360    Collective on Mat
3361 
3362    Input Parameters:
3363 +  mat - the matrix
3364 .  perm - row and column permutation
3365 .  fill - levels of fill
3366 -  f - expected fill as ratio of original fill
3367 
3368    Output Parameter:
3369 .  fact - the factored matrix
3370 
3371    Notes:
3372    Currently only no-fill factorization is supported.
3373 
3374    Most users should employ the simplified SLES interface for linear solvers
3375    instead of working directly with matrix algebra routines such as this.
3376    See, e.g., SLESCreate().
3377 
3378    Level: developer
3379 
3380   Concepts: matrices^symbolic incomplete Cholesky factorization
3381   Concepts: matrices^factorization
3382   Concepts: Cholsky^symbolic factorization
3383 
3384 .seealso: MatCholeskyFactorNumeric(), MatCholeskyFactor()
3385 @*/
3386 int MatIncompleteCholeskyFactorSymbolic(Mat mat,IS perm,PetscReal f,int fill,Mat *fact)
3387 {
3388   int ierr;
3389 
3390   PetscFunctionBegin;
3391   PetscValidHeaderSpecific(mat,MAT_COOKIE);
3392   PetscValidType(mat);
3393   MatPreallocated(mat);
3394   PetscValidPointer(fact);
3395   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
3396   if (fill < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Fill negative %d",fill);
3397   if (f < 1.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Expected fill less than 1.0 %g",f);
3398   if (!mat->ops->incompletecholeskyfactorsymbolic) SETERRQ1(PETSC_ERR_SUP,"Matrix type %s  symbolic ICC",mat->type_name);
3399   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
3400 
3401   ierr = PetscLogEventBegin(MAT_IncompleteCholeskyFactorSymbolic,mat,perm,0,0);CHKERRQ(ierr);
3402   ierr = (*mat->ops->incompletecholeskyfactorsymbolic)(mat,perm,f,fill,fact);CHKERRQ(ierr);
3403   ierr = PetscLogEventEnd(MAT_IncompleteCholeskyFactorSymbolic,mat,perm,0,0);CHKERRQ(ierr);
3404   PetscFunctionReturn(0);
3405 }
3406 
3407 #undef __FUNC__
3408 #define __FUNC__ "MatGetArray"
3409 /*@C
3410    MatGetArray - Returns a pointer to the element values in the matrix.
3411    The result of this routine is dependent on the underlying matrix data
3412    structure, and may not even work for certain matrix types.  You MUST
3413    call MatRestoreArray() when you no longer need to access the array.
3414 
3415    Not Collective
3416 
3417    Input Parameter:
3418 .  mat - the matrix
3419 
3420    Output Parameter:
3421 .  v - the location of the values
3422 
3423 
3424    Fortran Note:
3425    This routine is used differently from Fortran, e.g.,
3426 .vb
3427         Mat         mat
3428         Scalar      mat_array(1)
3429         PetscOffset i_mat
3430         int         ierr
3431         call MatGetArray(mat,mat_array,i_mat,ierr)
3432 
3433   C  Access first local entry in matrix; note that array is
3434   C  treated as one dimensional
3435         value = mat_array(i_mat + 1)
3436 
3437         [... other code ...]
3438         call MatRestoreArray(mat,mat_array,i_mat,ierr)
3439 .ve
3440 
3441    See the Fortran chapter of the users manual and
3442    petsc/src/mat/examples/tests for details.
3443 
3444    Level: advanced
3445 
3446    Concepts: matrices^access array
3447 
3448 .seealso: MatRestoreArray(), MatGetArrayF90()
3449 @*/
3450 int MatGetArray(Mat mat,Scalar **v)
3451 {
3452   int ierr;
3453 
3454   PetscFunctionBegin;
3455   PetscValidHeaderSpecific(mat,MAT_COOKIE);
3456   PetscValidType(mat);
3457   MatPreallocated(mat);
3458   PetscValidPointer(v);
3459   if (!mat->ops->getarray) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name);
3460   ierr = (*mat->ops->getarray)(mat,v);CHKERRQ(ierr);
3461   PetscFunctionReturn(0);
3462 }
3463 
3464 #undef __FUNC__
3465 #define __FUNC__ "MatRestoreArray"
3466 /*@C
3467    MatRestoreArray - Restores the matrix after MatGetArray() has been called.
3468 
3469    Not Collective
3470 
3471    Input Parameter:
3472 +  mat - the matrix
3473 -  v - the location of the values
3474 
3475    Fortran Note:
3476    This routine is used differently from Fortran, e.g.,
3477 .vb
3478         Mat         mat
3479         Scalar      mat_array(1)
3480         PetscOffset i_mat
3481         int         ierr
3482         call MatGetArray(mat,mat_array,i_mat,ierr)
3483 
3484   C  Access first local entry in matrix; note that array is
3485   C  treated as one dimensional
3486         value = mat_array(i_mat + 1)
3487 
3488         [... other code ...]
3489         call MatRestoreArray(mat,mat_array,i_mat,ierr)
3490 .ve
3491 
3492    See the Fortran chapter of the users manual and
3493    petsc/src/mat/examples/tests for details
3494 
3495    Level: advanced
3496 
3497 .seealso: MatGetArray(), MatRestoreArrayF90()
3498 @*/
3499 int MatRestoreArray(Mat mat,Scalar **v)
3500 {
3501   int ierr;
3502 
3503   PetscFunctionBegin;
3504   PetscValidHeaderSpecific(mat,MAT_COOKIE);
3505   PetscValidType(mat);
3506   MatPreallocated(mat);
3507   PetscValidPointer(v);
3508   if (!mat->ops->restorearray) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name);
3509   ierr = (*mat->ops->restorearray)(mat,v);CHKERRQ(ierr);
3510   PetscFunctionReturn(0);
3511 }
3512 
3513 #undef __FUNC__
3514 #define __FUNC__ "MatGetSubMatrices"
3515 /*@C
3516    MatGetSubMatrices - Extracts several submatrices from a matrix. If submat
3517    points to an array of valid matrices, they may be reused to store the new
3518    submatrices.
3519 
3520    Collective on Mat
3521 
3522    Input Parameters:
3523 +  mat - the matrix
3524 .  n   - the number of submatrixes to be extracted (on this processor, may be zero)
3525 .  irow, icol - index sets of rows and columns to extract
3526 -  scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
3527 
3528    Output Parameter:
3529 .  submat - the array of submatrices
3530 
3531    Notes:
3532    MatGetSubMatrices() can extract only sequential submatrices
3533    (from both sequential and parallel matrices). Use MatGetSubMatrix()
3534    to extract a parallel submatrix.
3535 
3536    When extracting submatrices from a parallel matrix, each processor can
3537    form a different submatrix by setting the rows and columns of its
3538    individual index sets according to the local submatrix desired.
3539 
3540    When finished using the submatrices, the user should destroy
3541    them with MatDestroySubMatrices().
3542 
3543    MAT_REUSE_MATRIX can only be used when the nonzero structure of the
3544    original matrix has not changed from that last call to MatGetSubMatrices().
3545 
3546    This routine creates the matrices submat; you should NOT create them before
3547    calling it.
3548 
3549    Fortran Note:
3550    The Fortran interface is slightly different from that given below; it
3551    requires one to pass in  as submat a Mat (integer) array of size at least m.
3552 
3553    Level: advanced
3554 
3555    Concepts: matrices^accessing submatrices
3556    Concepts: submatrices
3557 
3558 .seealso: MatDestroyMatrices(), MatGetSubMatrix(), MatGetRow(), MatGetDiagonal()
3559 @*/
3560 int MatGetSubMatrices(Mat mat,int n,IS *irow,IS *icol,MatReuse scall,Mat **submat)
3561 {
3562   int ierr;
3563 
3564   PetscFunctionBegin;
3565   PetscValidHeaderSpecific(mat,MAT_COOKIE);
3566   PetscValidType(mat);
3567   MatPreallocated(mat);
3568   if (!mat->ops->getsubmatrices) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name);
3569   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
3570 
3571   ierr = PetscLogEventBegin(MAT_GetSubMatrices,mat,0,0,0);CHKERRQ(ierr);
3572   ierr = (*mat->ops->getsubmatrices)(mat,n,irow,icol,scall,submat);CHKERRQ(ierr);
3573   ierr = PetscLogEventEnd(MAT_GetSubMatrices,mat,0,0,0);CHKERRQ(ierr);
3574 
3575   PetscFunctionReturn(0);
3576 }
3577 
3578 #undef __FUNC__
3579 #define __FUNC__ "MatDestroyMatrices"
3580 /*@C
3581    MatDestroyMatrices - Destroys a set of matrices obtained with MatGetSubMatrices().
3582 
3583    Collective on Mat
3584 
3585    Input Parameters:
3586 +  n - the number of local matrices
3587 -  mat - the matrices
3588 
3589    Level: advanced
3590 
3591 .seealso: MatGetSubMatrices()
3592 @*/
3593 int MatDestroyMatrices(int n,Mat **mat)
3594 {
3595   int ierr,i;
3596 
3597   PetscFunctionBegin;
3598   if (n < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Trying to destroy negative number of matrices %d",n);
3599   PetscValidPointer(mat);
3600   for (i=0; i<n; i++) {
3601     ierr = MatDestroy((*mat)[i]);CHKERRQ(ierr);
3602   }
3603   /* memory is allocated even if n = 0 */
3604   ierr = PetscFree(*mat);CHKERRQ(ierr);
3605   PetscFunctionReturn(0);
3606 }
3607 
3608 #undef __FUNC__
3609 #define __FUNC__ "MatIncreaseOverlap"
3610 /*@
3611    MatIncreaseOverlap - Given a set of submatrices indicated by index sets,
3612    replaces the index sets by larger ones that represent submatrices with
3613    additional overlap.
3614 
3615    Collective on Mat
3616 
3617    Input Parameters:
3618 +  mat - the matrix
3619 .  n   - the number of index sets
3620 .  is  - the array of pointers to index sets
3621 -  ov  - the additional overlap requested
3622 
3623    Level: developer
3624 
3625    Concepts: overlap
3626    Concepts: ASM^computing overlap
3627 
3628 .seealso: MatGetSubMatrices()
3629 @*/
3630 int MatIncreaseOverlap(Mat mat,int n,IS *is,int ov)
3631 {
3632   int ierr;
3633 
3634   PetscFunctionBegin;
3635   PetscValidHeaderSpecific(mat,MAT_COOKIE);
3636   PetscValidType(mat);
3637   MatPreallocated(mat);
3638   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
3639   if (mat->factor)     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
3640 
3641   if (!ov) PetscFunctionReturn(0);
3642   if (!mat->ops->increaseoverlap) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name);
3643   ierr = PetscLogEventBegin(MAT_IncreaseOverlap,mat,0,0,0);CHKERRQ(ierr);
3644   ierr = (*mat->ops->increaseoverlap)(mat,n,is,ov);CHKERRQ(ierr);
3645   ierr = PetscLogEventEnd(MAT_IncreaseOverlap,mat,0,0,0);CHKERRQ(ierr);
3646   PetscFunctionReturn(0);
3647 }
3648 
3649 #undef __FUNC__
3650 #define __FUNC__ "MatPrintHelp"
3651 /*@
3652    MatPrintHelp - Prints all the options for the matrix.
3653 
3654    Collective on Mat
3655 
3656    Input Parameter:
3657 .  mat - the matrix
3658 
3659    Options Database Keys:
3660 +  -help - Prints matrix options
3661 -  -h - Prints matrix options
3662 
3663    Level: developer
3664 
3665 .seealso: MatCreate(), MatCreateXXX()
3666 @*/
3667 int MatPrintHelp(Mat mat)
3668 {
3669   static PetscTruth called = PETSC_FALSE;
3670   int               ierr;
3671   MPI_Comm          comm;
3672 
3673   PetscFunctionBegin;
3674   PetscValidHeaderSpecific(mat,MAT_COOKIE);
3675   PetscValidType(mat);
3676   MatPreallocated(mat);
3677 
3678   comm = mat->comm;
3679   if (!called) {
3680     ierr = (*PetscHelpPrintf)(comm,"General matrix options:\n");CHKERRQ(ierr);
3681     ierr = (*PetscHelpPrintf)(comm,"  -mat_view_info: view basic matrix info during MatAssemblyEnd()\n");CHKERRQ(ierr);
3682     ierr = (*PetscHelpPrintf)(comm,"  -mat_view_info_detailed: view detailed matrix info during MatAssemblyEnd()\n");CHKERRQ(ierr);
3683     ierr = (*PetscHelpPrintf)(comm,"  -mat_view_draw: draw nonzero matrix structure during MatAssemblyEnd()\n");CHKERRQ(ierr);
3684     ierr = (*PetscHelpPrintf)(comm,"      -draw_pause <sec>: set seconds of display pause\n");CHKERRQ(ierr);
3685     ierr = (*PetscHelpPrintf)(comm,"      -display <name>: set alternate display\n");CHKERRQ(ierr);
3686     called = PETSC_TRUE;
3687   }
3688   if (mat->ops->printhelp) {
3689     ierr = (*mat->ops->printhelp)(mat);CHKERRQ(ierr);
3690   }
3691   PetscFunctionReturn(0);
3692 }
3693 
3694 #undef __FUNC__
3695 #define __FUNC__ "MatGetBlockSize"
3696 /*@
3697    MatGetBlockSize - Returns the matrix block size; useful especially for the
3698    block row and block diagonal formats.
3699 
3700    Not Collective
3701 
3702    Input Parameter:
3703 .  mat - the matrix
3704 
3705    Output Parameter:
3706 .  bs - block size
3707 
3708    Notes:
3709    Block diagonal formats are MATSEQBDIAG, MATMPIBDIAG.
3710    Block row formats are MATSEQBAIJ, MATMPIBAIJ
3711 
3712    Level: intermediate
3713 
3714    Concepts: matrices^block size
3715 
3716 .seealso: MatCreateSeqBAIJ(), MatCreateMPIBAIJ(), MatCreateSeqBDiag(), MatCreateMPIBDiag()
3717 @*/
3718 int MatGetBlockSize(Mat mat,int *bs)
3719 {
3720   int ierr;
3721 
3722   PetscFunctionBegin;
3723   PetscValidHeaderSpecific(mat,MAT_COOKIE);
3724   PetscValidType(mat);
3725   MatPreallocated(mat);
3726   PetscValidIntPointer(bs);
3727   if (!mat->ops->getblocksize) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name);
3728   ierr = (*mat->ops->getblocksize)(mat,bs);CHKERRQ(ierr);
3729   PetscFunctionReturn(0);
3730 }
3731 
3732 #undef __FUNC__
3733 #define __FUNC__ "MatGetRowIJ"
3734 /*@C
3735     MatGetRowIJ - Returns the compressed row storage i and j indices for sequential matrices.
3736 
3737    Collective on Mat
3738 
3739     Input Parameters:
3740 +   mat - the matrix
3741 .   shift -  0 or 1 indicating we want the indices starting at 0 or 1
3742 -   symmetric - PETSC_TRUE or PETSC_FALSE indicating the matrix data structure should be
3743                 symmetrized
3744 
3745     Output Parameters:
3746 +   n - number of rows in the (possibly compressed) matrix
3747 .   ia - the row pointers
3748 .   ja - the column indices
3749 -   done - PETSC_TRUE or PETSC_FALSE, indicating whether the values have been returned
3750 
3751     Level: developer
3752 
3753 .seealso: MatGetColumnIJ(), MatRestoreRowIJ()
3754 @*/
3755 int MatGetRowIJ(Mat mat,int shift,PetscTruth symmetric,int *n,int **ia,int** ja,PetscTruth *done)
3756 {
3757   int ierr;
3758 
3759   PetscFunctionBegin;
3760   PetscValidHeaderSpecific(mat,MAT_COOKIE);
3761   PetscValidType(mat);
3762   MatPreallocated(mat);
3763   if (ia) PetscValidIntPointer(ia);
3764   if (ja) PetscValidIntPointer(ja);
3765   PetscValidIntPointer(done);
3766   if (!mat->ops->getrowij) *done = PETSC_FALSE;
3767   else {
3768     *done = PETSC_TRUE;
3769     ierr  = (*mat->ops->getrowij)(mat,shift,symmetric,n,ia,ja,done);CHKERRQ(ierr);
3770   }
3771   PetscFunctionReturn(0);
3772 }
3773 
3774 #undef __FUNC__
3775 #define __FUNC__ "MatGetColumnIJ"
3776 /*@C
3777     MatGetColumnIJ - Returns the compressed column storage i and j indices for sequential matrices.
3778 
3779     Collective on Mat
3780 
3781     Input Parameters:
3782 +   mat - the matrix
3783 .   shift - 1 or zero indicating we want the indices starting at 0 or 1
3784 -   symmetric - PETSC_TRUE or PETSC_FALSE indicating the matrix data structure should be
3785                 symmetrized
3786 
3787     Output Parameters:
3788 +   n - number of columns in the (possibly compressed) matrix
3789 .   ia - the column pointers
3790 .   ja - the row indices
3791 -   done - PETSC_TRUE or PETSC_FALSE, indicating whether the values have been returned
3792 
3793     Level: developer
3794 
3795 .seealso: MatGetRowIJ(), MatRestoreColumnIJ()
3796 @*/
3797 int MatGetColumnIJ(Mat mat,int shift,PetscTruth symmetric,int *n,int **ia,int** ja,PetscTruth *done)
3798 {
3799   int ierr;
3800 
3801   PetscFunctionBegin;
3802   PetscValidHeaderSpecific(mat,MAT_COOKIE);
3803   PetscValidType(mat);
3804   MatPreallocated(mat);
3805   if (ia) PetscValidIntPointer(ia);
3806   if (ja) PetscValidIntPointer(ja);
3807   PetscValidIntPointer(done);
3808 
3809   if (!mat->ops->getcolumnij) *done = PETSC_FALSE;
3810   else {
3811     *done = PETSC_TRUE;
3812     ierr  = (*mat->ops->getcolumnij)(mat,shift,symmetric,n,ia,ja,done);CHKERRQ(ierr);
3813   }
3814   PetscFunctionReturn(0);
3815 }
3816 
3817 #undef __FUNC__
3818 #define __FUNC__ "MatRestoreRowIJ"
3819 /*@C
3820     MatRestoreRowIJ - Call after you are completed with the ia,ja indices obtained with
3821     MatGetRowIJ().
3822 
3823     Collective on Mat
3824 
3825     Input Parameters:
3826 +   mat - the matrix
3827 .   shift - 1 or zero indicating we want the indices starting at 0 or 1
3828 -   symmetric - PETSC_TRUE or PETSC_FALSE indicating the matrix data structure should be
3829                 symmetrized
3830 
3831     Output Parameters:
3832 +   n - size of (possibly compressed) matrix
3833 .   ia - the row pointers
3834 .   ja - the column indices
3835 -   done - PETSC_TRUE or PETSC_FALSE indicated that the values have been returned
3836 
3837     Level: developer
3838 
3839 .seealso: MatGetRowIJ(), MatRestoreColumnIJ()
3840 @*/
3841 int MatRestoreRowIJ(Mat mat,int shift,PetscTruth symmetric,int *n,int **ia,int** ja,PetscTruth *done)
3842 {
3843   int ierr;
3844 
3845   PetscFunctionBegin;
3846   PetscValidHeaderSpecific(mat,MAT_COOKIE);
3847   PetscValidType(mat);
3848   MatPreallocated(mat);
3849   if (ia) PetscValidIntPointer(ia);
3850   if (ja) PetscValidIntPointer(ja);
3851   PetscValidIntPointer(done);
3852 
3853   if (!mat->ops->restorerowij) *done = PETSC_FALSE;
3854   else {
3855     *done = PETSC_TRUE;
3856     ierr  = (*mat->ops->restorerowij)(mat,shift,symmetric,n,ia,ja,done);CHKERRQ(ierr);
3857   }
3858   PetscFunctionReturn(0);
3859 }
3860 
3861 #undef __FUNC__
3862 #define __FUNC__ "MatRestoreColumnIJ"
3863 /*@C
3864     MatRestoreColumnIJ - Call after you are completed with the ia,ja indices obtained with
3865     MatGetColumnIJ().
3866 
3867     Collective on Mat
3868 
3869     Input Parameters:
3870 +   mat - the matrix
3871 .   shift - 1 or zero indicating we want the indices starting at 0 or 1
3872 -   symmetric - PETSC_TRUE or PETSC_FALSE indicating the matrix data structure should be
3873                 symmetrized
3874 
3875     Output Parameters:
3876 +   n - size of (possibly compressed) matrix
3877 .   ia - the column pointers
3878 .   ja - the row indices
3879 -   done - PETSC_TRUE or PETSC_FALSE indicated that the values have been returned
3880 
3881     Level: developer
3882 
3883 .seealso: MatGetColumnIJ(), MatRestoreRowIJ()
3884 @*/
3885 int MatRestoreColumnIJ(Mat mat,int shift,PetscTruth symmetric,int *n,int **ia,int** ja,PetscTruth *done)
3886 {
3887   int ierr;
3888 
3889   PetscFunctionBegin;
3890   PetscValidHeaderSpecific(mat,MAT_COOKIE);
3891   PetscValidType(mat);
3892   MatPreallocated(mat);
3893   if (ia) PetscValidIntPointer(ia);
3894   if (ja) PetscValidIntPointer(ja);
3895   PetscValidIntPointer(done);
3896 
3897   if (!mat->ops->restorecolumnij) *done = PETSC_FALSE;
3898   else {
3899     *done = PETSC_TRUE;
3900     ierr  = (*mat->ops->restorecolumnij)(mat,shift,symmetric,n,ia,ja,done);CHKERRQ(ierr);
3901   }
3902   PetscFunctionReturn(0);
3903 }
3904 
3905 #undef __FUNC__
3906 #define __FUNC__ "MatColoringPatch"
3907 /*@C
3908     MatColoringPatch -Used inside matrix coloring routines that
3909     use MatGetRowIJ() and/or MatGetColumnIJ().
3910 
3911     Collective on Mat
3912 
3913     Input Parameters:
3914 +   mat - the matrix
3915 .   n   - number of colors
3916 -   colorarray - array indicating color for each column
3917 
3918     Output Parameters:
3919 .   iscoloring - coloring generated using colorarray information
3920 
3921     Level: developer
3922 
3923 .seealso: MatGetRowIJ(), MatGetColumnIJ()
3924 
3925 @*/
3926 int MatColoringPatch(Mat mat,int n,int *colorarray,ISColoring *iscoloring)
3927 {
3928   int ierr;
3929 
3930   PetscFunctionBegin;
3931   PetscValidHeaderSpecific(mat,MAT_COOKIE);
3932   PetscValidType(mat);
3933   MatPreallocated(mat);
3934   PetscValidIntPointer(colorarray);
3935 
3936   if (!mat->ops->coloringpatch){
3937     SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name);
3938   } else {
3939     ierr  = (*mat->ops->coloringpatch)(mat,n,colorarray,iscoloring);CHKERRQ(ierr);
3940   }
3941   PetscFunctionReturn(0);
3942 }
3943 
3944 
3945 #undef __FUNC__
3946 #define __FUNC__ "MatSetUnfactored"
3947 /*@
3948    MatSetUnfactored - Resets a factored matrix to be treated as unfactored.
3949 
3950    Collective on Mat
3951 
3952    Input Parameter:
3953 .  mat - the factored matrix to be reset
3954 
3955    Notes:
3956    This routine should be used only with factored matrices formed by in-place
3957    factorization via ILU(0) (or by in-place LU factorization for the MATSEQDENSE
3958    format).  This option can save memory, for example, when solving nonlinear
3959    systems with a matrix-free Newton-Krylov method and a matrix-based, in-place
3960    ILU(0) preconditioner.
3961 
3962    Note that one can specify in-place ILU(0) factorization by calling
3963 .vb
3964      PCType(pc,PCILU);
3965      PCILUSeUseInPlace(pc);
3966 .ve
3967    or by using the options -pc_type ilu -pc_ilu_in_place
3968 
3969    In-place factorization ILU(0) can also be used as a local
3970    solver for the blocks within the block Jacobi or additive Schwarz
3971    methods (runtime option: -sub_pc_ilu_in_place).  See the discussion
3972    of these preconditioners in the users manual for details on setting
3973    local solver options.
3974 
3975    Most users should employ the simplified SLES interface for linear solvers
3976    instead of working directly with matrix algebra routines such as this.
3977    See, e.g., SLESCreate().
3978 
3979    Level: developer
3980 
3981 .seealso: PCILUSetUseInPlace(), PCLUSetUseInPlace()
3982 
3983    Concepts: matrices^unfactored
3984 
3985 @*/
3986 int MatSetUnfactored(Mat mat)
3987 {
3988   int ierr;
3989 
3990   PetscFunctionBegin;
3991   PetscValidHeaderSpecific(mat,MAT_COOKIE);
3992   PetscValidType(mat);
3993   MatPreallocated(mat);
3994   mat->factor = 0;
3995   if (!mat->ops->setunfactored) PetscFunctionReturn(0);
3996   ierr = (*mat->ops->setunfactored)(mat);CHKERRQ(ierr);
3997   PetscFunctionReturn(0);
3998 }
3999 
4000 /*MC
4001     MatGetArrayF90 - Accesses a matrix array from Fortran90.
4002 
4003     Synopsis:
4004     MatGetArrayF90(Mat x,{Scalar, pointer :: xx_v(:)},integer ierr)
4005 
4006     Not collective
4007 
4008     Input Parameter:
4009 .   x - matrix
4010 
4011     Output Parameters:
4012 +   xx_v - the Fortran90 pointer to the array
4013 -   ierr - error code
4014 
4015     Example of Usage:
4016 .vb
4017       Scalar, pointer xx_v(:)
4018       ....
4019       call MatGetArrayF90(x,xx_v,ierr)
4020       a = xx_v(3)
4021       call MatRestoreArrayF90(x,xx_v,ierr)
4022 .ve
4023 
4024     Notes:
4025     Not yet supported for all F90 compilers
4026 
4027     Level: advanced
4028 
4029 .seealso:  MatRestoreArrayF90(), MatGetArray(), MatRestoreArray()
4030 
4031     Concepts: matrices^accessing array
4032 
4033 M*/
4034 
4035 /*MC
4036     MatRestoreArrayF90 - Restores a matrix array that has been
4037     accessed with MatGetArrayF90().
4038 
4039     Synopsis:
4040     MatRestoreArrayF90(Mat x,{Scalar, pointer :: xx_v(:)},integer ierr)
4041 
4042     Not collective
4043 
4044     Input Parameters:
4045 +   x - matrix
4046 -   xx_v - the Fortran90 pointer to the array
4047 
4048     Output Parameter:
4049 .   ierr - error code
4050 
4051     Example of Usage:
4052 .vb
4053        Scalar, pointer xx_v(:)
4054        ....
4055        call MatGetArrayF90(x,xx_v,ierr)
4056        a = xx_v(3)
4057        call MatRestoreArrayF90(x,xx_v,ierr)
4058 .ve
4059 
4060     Notes:
4061     Not yet supported for all F90 compilers
4062 
4063     Level: advanced
4064 
4065 .seealso:  MatGetArrayF90(), MatGetArray(), MatRestoreArray()
4066 
4067 M*/
4068 
4069 
4070 #undef __FUNC__
4071 #define __FUNC__ "MatGetSubMatrix"
4072 /*@
4073     MatGetSubMatrix - Gets a single submatrix on the same number of processors
4074                       as the original matrix.
4075 
4076     Collective on Mat
4077 
4078     Input Parameters:
4079 +   mat - the original matrix
4080 .   isrow - rows this processor should obtain
4081 .   iscol - columns for all processors you wish to keep
4082 .   csize - number of columns "local" to this processor (does nothing for sequential
4083             matrices). This should match the result from VecGetLocalSize(x,...) if you
4084             plan to use the matrix in a A*x; alternatively, you can use PETSC_DECIDE
4085 -   cll - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
4086 
4087     Output Parameter:
4088 .   newmat - the new submatrix, of the same type as the old
4089 
4090     Level: advanced
4091 
4092     Notes: the iscol argument MUST be the same on each processor.
4093 
4094       The first time this is called you should use a cll of MAT_INITIAL_MATRIX,
4095    the MatGetSubMatrix() routine will create the newmat for you. Any additional calls
4096    to this routine with a mat of the same nonzero structure will reuse the matrix
4097    generated the first time.
4098 
4099     Concepts: matrices^submatrices
4100 
4101 .seealso: MatGetSubMatrices()
4102 @*/
4103 int MatGetSubMatrix(Mat mat,IS isrow,IS iscol,int csize,MatReuse cll,Mat *newmat)
4104 {
4105   int     ierr, size;
4106   Mat     *local;
4107 
4108   PetscFunctionBegin;
4109   PetscValidType(mat);
4110   MatPreallocated(mat);
4111   ierr = MPI_Comm_size(mat->comm,&size);CHKERRQ(ierr);
4112 
4113   /* if original matrix is on just one processor then use submatrix generated */
4114   if (!mat->ops->getsubmatrix && size == 1 && cll == MAT_REUSE_MATRIX) {
4115     ierr = MatGetSubMatrices(mat,1,&isrow,&iscol,MAT_REUSE_MATRIX,&newmat);CHKERRQ(ierr);
4116     PetscFunctionReturn(0);
4117   } else if (!mat->ops->getsubmatrix && size == 1) {
4118     ierr = MatGetSubMatrices(mat,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&local);CHKERRQ(ierr);
4119     *newmat = *local;
4120     ierr = PetscFree(local);CHKERRQ(ierr);
4121     PetscFunctionReturn(0);
4122   }
4123 
4124   if (!mat->ops->getsubmatrix) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name);
4125   ierr = (*mat->ops->getsubmatrix)(mat,isrow,iscol,csize,cll,newmat);CHKERRQ(ierr);
4126   PetscFunctionReturn(0);
4127 }
4128 
4129 #undef __FUNC__
4130 #define __FUNC__ "MatGetMaps"
4131 /*@C
4132    MatGetMaps - Returns the maps associated with the matrix.
4133 
4134    Not Collective
4135 
4136    Input Parameter:
4137 .  mat - the matrix
4138 
4139    Output Parameters:
4140 +  rmap - the row (right) map
4141 -  cmap - the column (left) map
4142 
4143    Level: developer
4144 
4145    Concepts: maps^getting from matrix
4146 
4147 @*/
4148 int MatGetMaps(Mat mat,Map *rmap,Map *cmap)
4149 {
4150   int ierr;
4151 
4152   PetscFunctionBegin;
4153   PetscValidHeaderSpecific(mat,MAT_COOKIE);
4154   PetscValidType(mat);
4155   MatPreallocated(mat);
4156   ierr = (*mat->ops->getmaps)(mat,rmap,cmap);CHKERRQ(ierr);
4157   PetscFunctionReturn(0);
4158 }
4159 
4160 /*
4161       Version that works for all PETSc matrices
4162 */
4163 #undef __FUNC__
4164 #define __FUNC__ "MatGetMaps_Petsc"
4165 int MatGetMaps_Petsc(Mat mat,Map *rmap,Map *cmap)
4166 {
4167   PetscFunctionBegin;
4168   if (rmap) *rmap = mat->rmap;
4169   if (cmap) *cmap = mat->cmap;
4170   PetscFunctionReturn(0);
4171 }
4172 
4173 #undef __FUNC__
4174 #define __FUNC__ "MatSetStashInitialSize"
4175 /*@
4176    MatSetStashInitialSize - sets the sizes of the matrix stash, that is
4177    used during the assembly process to store values that belong to
4178    other processors.
4179 
4180    Not Collective
4181 
4182    Input Parameters:
4183 +  mat   - the matrix
4184 .  size  - the initial size of the stash.
4185 -  bsize - the initial size of the block-stash(if used).
4186 
4187    Options Database Keys:
4188 +   -matstash_initial_size <size> or <size0,size1,...sizep-1>
4189 -   -matstash_block_initial_size <bsize>  or <bsize0,bsize1,...bsizep-1>
4190 
4191    Level: intermediate
4192 
4193    Notes:
4194      The block-stash is used for values set with VecSetValuesBlocked() while
4195      the stash is used for values set with VecSetValues()
4196 
4197      Run with the option -log_info and look for output of the form
4198      MatAssemblyBegin_MPIXXX:Stash has MM entries, uses nn mallocs.
4199      to determine the appropriate value, MM, to use for size and
4200      MatAssemblyBegin_MPIXXX:Block-Stash has BMM entries, uses nn mallocs.
4201      to determine the value, BMM to use for bsize
4202 
4203    Concepts: stash^setting matrix size
4204    Concepts: matrices^stash
4205 
4206 @*/
4207 int MatSetStashInitialSize(Mat mat,int size, int bsize)
4208 {
4209   int ierr;
4210 
4211   PetscFunctionBegin;
4212   PetscValidHeaderSpecific(mat,MAT_COOKIE);
4213   PetscValidType(mat);
4214   MatPreallocated(mat);
4215   ierr = MatStashSetInitialSize_Private(&mat->stash,size);CHKERRQ(ierr);
4216   ierr = MatStashSetInitialSize_Private(&mat->bstash,bsize);CHKERRQ(ierr);
4217   PetscFunctionReturn(0);
4218 }
4219 
4220 #undef __FUNC__
4221 #define __FUNC__ "MatInterpolateAdd"
4222 /*@
4223    MatInterpolateAdd - w = y + A*x or A'*x depending on the shape of
4224      the matrix
4225 
4226    Collective on Mat
4227 
4228    Input Parameters:
4229 +  mat   - the matrix
4230 .  x,y - the vectors
4231 -  w - where the result is stored
4232 
4233    Level: intermediate
4234 
4235    Notes:
4236     w may be the same vector as y.
4237 
4238     This allows one to use either the restriction or interpolation (its transpose)
4239     matrix to do the interpolation
4240 
4241     Concepts: interpolation
4242 
4243 .seealso: MatMultAdd(), MatMultTransposeAdd(), MatRestrict()
4244 
4245 @*/
4246 int MatInterpolateAdd(Mat A,Vec x,Vec y,Vec w)
4247 {
4248   int M,N,ierr;
4249 
4250   PetscFunctionBegin;
4251   PetscValidType(A);
4252   MatPreallocated(A);
4253   ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr);
4254   if (N > M) {
4255     ierr = MatMultTransposeAdd(A,x,y,w);CHKERRQ(ierr);
4256   } else {
4257     ierr = MatMultAdd(A,x,y,w);CHKERRQ(ierr);
4258   }
4259   PetscFunctionReturn(0);
4260 }
4261 
4262 #undef __FUNC__
4263 #define __FUNC__ "MatInterpolate"
4264 /*@
4265    MatInterpolate - y = A*x or A'*x depending on the shape of
4266      the matrix
4267 
4268    Collective on Mat
4269 
4270    Input Parameters:
4271 +  mat   - the matrix
4272 -  x,y - the vectors
4273 
4274    Level: intermediate
4275 
4276    Notes:
4277     This allows one to use either the restriction or interpolation (its transpose)
4278     matrix to do the interpolation
4279 
4280    Concepts: matrices^interpolation
4281 
4282 .seealso: MatMultAdd(), MatMultTransposeAdd(), MatRestrict()
4283 
4284 @*/
4285 int MatInterpolate(Mat A,Vec x,Vec y)
4286 {
4287   int M,N,ierr;
4288 
4289   PetscFunctionBegin;
4290   PetscValidType(A);
4291   MatPreallocated(A);
4292   ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr);
4293   if (N > M) {
4294     ierr = MatMultTranspose(A,x,y);CHKERRQ(ierr);
4295   } else {
4296     ierr = MatMult(A,x,y);CHKERRQ(ierr);
4297   }
4298   PetscFunctionReturn(0);
4299 }
4300 
4301 #undef __FUNC__
4302 #define __FUNC__ "MatRestrict"
4303 /*@
4304    MatRestrict - y = A*x or A'*x
4305 
4306    Collective on Mat
4307 
4308    Input Parameters:
4309 +  mat   - the matrix
4310 -  x,y - the vectors
4311 
4312    Level: intermediate
4313 
4314    Notes:
4315     This allows one to use either the restriction or interpolation (its transpose)
4316     matrix to do the restriction
4317 
4318    Concepts: matrices^restriction
4319 
4320 .seealso: MatMultAdd(), MatMultTransposeAdd(), MatInterpolate()
4321 
4322 @*/
4323 int MatRestrict(Mat A,Vec x,Vec y)
4324 {
4325   int M,N,ierr;
4326 
4327   PetscFunctionBegin;
4328   PetscValidType(A);
4329   MatPreallocated(A);
4330   ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr);
4331   if (N > M) {
4332     ierr = MatMult(A,x,y);CHKERRQ(ierr);
4333   } else {
4334     ierr = MatMultTranspose(A,x,y);CHKERRQ(ierr);
4335   }
4336   PetscFunctionReturn(0);
4337 }
4338 
4339 #undef __FUNC__
4340 #define __FUNC__ "MatNullSpaceAttach"
4341 /*@C
4342    MatNullSpaceAttach - attaches a null space to a matrix.
4343         This null space will be removed from the resulting vector whenever
4344         MatMult() is called
4345 
4346    Collective on Mat
4347 
4348    Input Parameters:
4349 +  mat - the matrix
4350 -  nullsp - the null space object
4351 
4352    Level: developer
4353 
4354    Notes:
4355       Overwrites any previous null space that may have been attached
4356 
4357    Concepts: null space^attaching to matrix
4358 
4359 .seealso: MatCreate(), MatNullSpaceCreate()
4360 @*/
4361 int MatNullSpaceAttach(Mat mat,MatNullSpace nullsp)
4362 {
4363   int ierr;
4364 
4365   PetscFunctionBegin;
4366   PetscValidHeaderSpecific(mat,MAT_COOKIE);
4367   PetscValidType(mat);
4368   MatPreallocated(mat);
4369   PetscValidHeaderSpecific(nullsp,MATNULLSPACE_COOKIE);
4370 
4371   if (mat->nullsp) {
4372     ierr = MatNullSpaceDestroy(mat->nullsp);CHKERRQ(ierr);
4373   }
4374   mat->nullsp = nullsp;
4375   ierr = PetscObjectReference((PetscObject)nullsp);CHKERRQ(ierr);
4376   PetscFunctionReturn(0);
4377 }
4378 
4379 #undef __FUNC__
4380 #define __FUNC__ "MatIncompleteCholeskyFactor"
4381 /*@
4382    MatIncompleteCholeskyFactor - Performs in-place incomplete Cholesky factorization of matrix.
4383 
4384    Collective on Mat
4385 
4386    Input Parameters:
4387 +  mat - the matrix
4388 .  row - row/column permutation
4389 .  fill - expected fill factor >= 1.0
4390 -  level - level of fill, for ICC(k)
4391 
4392    Notes:
4393    Probably really in-place only when level of fill is zero, otherwise allocates
4394    new space to store factored matrix and deletes previous memory.
4395 
4396    Most users should employ the simplified SLES interface for linear solvers
4397    instead of working directly with matrix algebra routines such as this.
4398    See, e.g., SLESCreate().
4399 
4400    Level: developer
4401 
4402    Concepts: matrices^incomplete Cholesky factorization
4403    Concepts: Cholesky factorization
4404 
4405 .seealso: MatIncompleteCholeskyFactorSymbolic(), MatLUFactorNumeric(), MatCholeskyFactor()
4406 @*/
4407 int MatIncompleteCholeskyFactor(Mat mat,IS row,PetscReal fill,int level)
4408 {
4409   int ierr;
4410 
4411   PetscFunctionBegin;
4412   PetscValidHeaderSpecific(mat,MAT_COOKIE);
4413   PetscValidType(mat);
4414   MatPreallocated(mat);
4415   if (mat->M != mat->N) SETERRQ(PETSC_ERR_ARG_WRONG,"matrix must be square");
4416   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
4417   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
4418   if (!mat->ops->incompletecholeskyfactor) SETERRQ1(PETSC_ERR_SUP,"Mat type %s",mat->type_name);
4419   ierr = (*mat->ops->incompletecholeskyfactor)(mat,row,fill,level);CHKERRQ(ierr);
4420   PetscFunctionReturn(0);
4421 }
4422