xref: /petsc/src/mat/impls/shell/shell.c (revision 63f3c55c12ae2f190c391f6df6c540efe018a44a)
1 
2 /*
3    This provides a simple shell for Fortran (and C programmers) to
4   create a very simple matrix class for use with KSP without coding
5   much of anything.
6 */
7 
8 #include <petsc/private/matimpl.h>        /*I "petscmat.h" I*/
9 #include <petsc/private/vecimpl.h>
10 
11 struct _MatShellOps {
12   /*  3 */ PetscErrorCode (*mult)(Mat,Vec,Vec);
13   /*  5 */ PetscErrorCode (*multtranspose)(Mat,Vec,Vec);
14   /* 17 */ PetscErrorCode (*getdiagonal)(Mat,Vec);
15   /* 43 */ PetscErrorCode (*copy)(Mat,Mat,MatStructure);
16   /* 60 */ PetscErrorCode (*destroy)(Mat);
17 };
18 
19 typedef struct {
20   struct _MatShellOps ops[1];
21 
22   PetscScalar vscale,vshift;
23   Vec         dshift;
24   Vec         left,right;
25   Vec         left_work,right_work;
26   Vec         left_add_work,right_add_work;
27   Mat         axpy;
28   PetscScalar axpy_vscale;
29   PetscBool   managescalingshifts;                   /* The user will manage the scaling and shifts for the MATSHELL, not the default */
30   /* support for ZeroRows/Columns operations */
31   IS          zrows;
32   IS          zcols;
33   Vec         zvals;
34   Vec         zvals_w;
35   VecScatter  zvals_sct_r;
36   VecScatter  zvals_sct_c;
37   void        *ctx;
38 } Mat_Shell;
39 
40 
41 /*
42      Store and scale values on zeroed rows
43      xx = [x_1, 0], 0 on zeroed columns
44 */
45 static PetscErrorCode MatShellPreZeroRight(Mat A,Vec x,Vec *xx)
46 {
47   Mat_Shell      *shell = (Mat_Shell*)A->data;
48   PetscErrorCode ierr;
49 
50   PetscFunctionBegin;
51   *xx = x;
52   if (shell->zrows) {
53     ierr = VecSet(shell->zvals_w,0.0);CHKERRQ(ierr);
54     ierr = VecScatterBegin(shell->zvals_sct_c,x,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
55     ierr = VecScatterEnd(shell->zvals_sct_c,x,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
56     ierr = VecPointwiseMult(shell->zvals_w,shell->zvals_w,shell->zvals);CHKERRQ(ierr);
57   }
58   if (shell->zcols) {
59     if (!shell->right_work) {
60       ierr = MatCreateVecs(A,&shell->right_work,NULL);CHKERRQ(ierr);
61     }
62     ierr = VecCopy(x,shell->right_work);CHKERRQ(ierr);
63     ierr = VecISSet(shell->right_work,shell->zcols,0.0);CHKERRQ(ierr);
64     *xx  = shell->right_work;
65   }
66   PetscFunctionReturn(0);
67 }
68 
69 /* Insert properly diagonally scaled values stored in MatShellPreZeroRight */
70 static PetscErrorCode MatShellPostZeroLeft(Mat A,Vec x)
71 {
72   Mat_Shell      *shell = (Mat_Shell*)A->data;
73   PetscErrorCode ierr;
74 
75   PetscFunctionBegin;
76   if (shell->zrows) {
77     ierr = VecScatterBegin(shell->zvals_sct_r,shell->zvals_w,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
78     ierr = VecScatterEnd(shell->zvals_sct_r,shell->zvals_w,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
79   }
80   PetscFunctionReturn(0);
81 }
82 
83 /*
84      Store and scale values on zeroed rows
85      xx = [x_1, 0], 0 on zeroed rows
86 */
87 static PetscErrorCode MatShellPreZeroLeft(Mat A,Vec x,Vec *xx)
88 {
89   Mat_Shell      *shell = (Mat_Shell*)A->data;
90   PetscErrorCode ierr;
91 
92   PetscFunctionBegin;
93   *xx = NULL;
94   if (!shell->zrows) {
95     *xx = x;
96   } else {
97     if (!shell->left_work) {
98       ierr = MatCreateVecs(A,NULL,&shell->left_work);CHKERRQ(ierr);
99     }
100     ierr = VecCopy(x,shell->left_work);CHKERRQ(ierr);
101     ierr = VecSet(shell->zvals_w,0.0);CHKERRQ(ierr);
102     ierr = VecScatterBegin(shell->zvals_sct_r,shell->zvals_w,shell->left_work,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
103     ierr = VecScatterEnd(shell->zvals_sct_r,shell->zvals_w,shell->left_work,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
104     ierr = VecScatterBegin(shell->zvals_sct_r,x,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
105     ierr = VecScatterEnd(shell->zvals_sct_r,x,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
106     ierr = VecPointwiseMult(shell->zvals_w,shell->zvals_w,shell->zvals);CHKERRQ(ierr);
107     *xx  = shell->left_work;
108   }
109   PetscFunctionReturn(0);
110 }
111 
112 /* Zero zero-columns contributions, sum contributions from properly scaled values stored in MatShellPreZeroLeft */
113 static PetscErrorCode MatShellPostZeroRight(Mat A,Vec x)
114 {
115   Mat_Shell      *shell = (Mat_Shell*)A->data;
116   PetscErrorCode ierr;
117 
118   PetscFunctionBegin;
119   if (shell->zcols) {
120     ierr = VecISSet(x,shell->zcols,0.0);CHKERRQ(ierr);
121   }
122   if (shell->zrows) {
123     ierr = VecScatterBegin(shell->zvals_sct_c,shell->zvals_w,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
124     ierr = VecScatterEnd(shell->zvals_sct_c,shell->zvals_w,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
125   }
126   PetscFunctionReturn(0);
127 }
128 
129 /*
130       xx = diag(left)*x
131 */
132 static PetscErrorCode MatShellPreScaleLeft(Mat A,Vec x,Vec *xx)
133 {
134   Mat_Shell      *shell = (Mat_Shell*)A->data;
135   PetscErrorCode ierr;
136 
137   PetscFunctionBegin;
138   *xx = NULL;
139   if (!shell->left) {
140     *xx = x;
141   } else {
142     if (!shell->left_work) {ierr = VecDuplicate(shell->left,&shell->left_work);CHKERRQ(ierr);}
143     ierr = VecPointwiseMult(shell->left_work,x,shell->left);CHKERRQ(ierr);
144     *xx  = shell->left_work;
145   }
146   PetscFunctionReturn(0);
147 }
148 
149 /*
150      xx = diag(right)*x
151 */
152 static PetscErrorCode MatShellPreScaleRight(Mat A,Vec x,Vec *xx)
153 {
154   Mat_Shell      *shell = (Mat_Shell*)A->data;
155   PetscErrorCode ierr;
156 
157   PetscFunctionBegin;
158   *xx = NULL;
159   if (!shell->right) {
160     *xx = x;
161   } else {
162     if (!shell->right_work) {ierr = VecDuplicate(shell->right,&shell->right_work);CHKERRQ(ierr);}
163     ierr = VecPointwiseMult(shell->right_work,x,shell->right);CHKERRQ(ierr);
164     *xx  = shell->right_work;
165   }
166   PetscFunctionReturn(0);
167 }
168 
169 /*
170     x = diag(left)*x
171 */
172 static PetscErrorCode MatShellPostScaleLeft(Mat A,Vec x)
173 {
174   Mat_Shell      *shell = (Mat_Shell*)A->data;
175   PetscErrorCode ierr;
176 
177   PetscFunctionBegin;
178   if (shell->left) {ierr = VecPointwiseMult(x,x,shell->left);CHKERRQ(ierr);}
179   PetscFunctionReturn(0);
180 }
181 
182 /*
183     x = diag(right)*x
184 */
185 static PetscErrorCode MatShellPostScaleRight(Mat A,Vec x)
186 {
187   Mat_Shell      *shell = (Mat_Shell*)A->data;
188   PetscErrorCode ierr;
189 
190   PetscFunctionBegin;
191   if (shell->right) {ierr = VecPointwiseMult(x,x,shell->right);CHKERRQ(ierr);}
192   PetscFunctionReturn(0);
193 }
194 
195 /*
196          Y = vscale*Y + diag(dshift)*X + vshift*X
197 
198          On input Y already contains A*x
199 */
200 static PetscErrorCode MatShellShiftAndScale(Mat A,Vec X,Vec Y)
201 {
202   Mat_Shell      *shell = (Mat_Shell*)A->data;
203   PetscErrorCode ierr;
204 
205   PetscFunctionBegin;
206   if (shell->dshift) {          /* get arrays because there is no VecPointwiseMultAdd() */
207     PetscInt          i,m;
208     const PetscScalar *x,*d;
209     PetscScalar       *y;
210     ierr = VecGetLocalSize(X,&m);CHKERRQ(ierr);
211     ierr = VecGetArrayRead(shell->dshift,&d);CHKERRQ(ierr);
212     ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
213     ierr = VecGetArray(Y,&y);CHKERRQ(ierr);
214     for (i=0; i<m; i++) y[i] = shell->vscale*y[i] + d[i]*x[i];
215     ierr = VecRestoreArrayRead(shell->dshift,&d);CHKERRQ(ierr);
216     ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
217     ierr = VecRestoreArray(Y,&y);CHKERRQ(ierr);
218   } else {
219     ierr = VecScale(Y,shell->vscale);CHKERRQ(ierr);
220   }
221   if (shell->vshift != 0.0) {ierr = VecAXPY(Y,shell->vshift,X);CHKERRQ(ierr);} /* if test is for non-square matrices */
222   PetscFunctionReturn(0);
223 }
224 
225 PetscErrorCode  MatShellGetContext_Shell(Mat mat,void *ctx)
226 {
227   PetscFunctionBegin;
228   *(void**)ctx = ((Mat_Shell*)(mat->data))->ctx;
229   PetscFunctionReturn(0);
230 }
231 
232 /*@
233     MatShellGetContext - Returns the user-provided context associated with a shell matrix.
234 
235     Not Collective
236 
237     Input Parameter:
238 .   mat - the matrix, should have been created with MatCreateShell()
239 
240     Output Parameter:
241 .   ctx - the user provided context
242 
243     Level: advanced
244 
245    Fortran Notes:
246     To use this from Fortran you must write a Fortran interface definition for this
247     function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.
248 
249 .seealso: MatCreateShell(), MatShellSetOperation(), MatShellSetContext()
250 @*/
251 PetscErrorCode  MatShellGetContext(Mat mat,void *ctx)
252 {
253   PetscErrorCode ierr;
254 
255   PetscFunctionBegin;
256   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
257   PetscValidPointer(ctx,2);
258   ierr = PetscUseMethod(mat,"MatShellGetContext_C",(Mat,void*),(mat,ctx));CHKERRQ(ierr);
259   PetscFunctionReturn(0);
260 }
261 
262 static PetscErrorCode MatZeroRowsColumns_Local_Shell(Mat mat,PetscInt nr,PetscInt rows[],PetscInt nc,PetscInt cols[],PetscScalar diag,PetscBool rc)
263 {
264   PetscErrorCode ierr;
265   Mat_Shell      *shell = (Mat_Shell*)mat->data;
266   Vec            x = NULL,b = NULL;
267   IS             is1, is2;
268   const PetscInt *ridxs;
269   PetscInt       *idxs,*gidxs;
270   PetscInt       cum,rst,cst,i;
271 
272   PetscFunctionBegin;
273   if (shell->axpy) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot diagonal scale MATSHELL after MatAXPY operation");
274   if (!shell->zvals) {
275     ierr = MatCreateVecs(mat,NULL,&shell->zvals);CHKERRQ(ierr);
276   }
277   if (!shell->zvals_w) {
278     ierr = VecDuplicate(shell->zvals,&shell->zvals_w);CHKERRQ(ierr);
279   }
280   ierr = MatGetOwnershipRange(mat,&rst,NULL);CHKERRQ(ierr);
281   ierr = MatGetOwnershipRangeColumn(mat,&cst,NULL);CHKERRQ(ierr);
282 
283   /* Expand/create index set of zeroed rows */
284   ierr = PetscMalloc1(nr,&idxs);CHKERRQ(ierr);
285   for (i = 0; i < nr; i++) idxs[i] = rows[i] + rst;
286   ierr = ISCreateGeneral(PETSC_COMM_SELF,nr,idxs,PETSC_OWN_POINTER,&is1);CHKERRQ(ierr);
287   ierr = ISSort(is1);CHKERRQ(ierr);
288   ierr = VecISSet(shell->zvals,is1,diag);CHKERRQ(ierr);
289   if (shell->zrows) {
290     ierr = ISSum(shell->zrows,is1,&is2);CHKERRQ(ierr);
291     ierr = ISDestroy(&shell->zrows);CHKERRQ(ierr);
292     ierr = ISDestroy(&is1);CHKERRQ(ierr);
293     shell->zrows = is2;
294   } else shell->zrows = is1;
295 
296   /* Create scatters for diagonal values communications */
297   ierr = VecScatterDestroy(&shell->zvals_sct_c);CHKERRQ(ierr);
298   ierr = VecScatterDestroy(&shell->zvals_sct_r);CHKERRQ(ierr);
299 
300   /* row scatter: from/to left vector */
301   ierr = MatCreateVecs(mat,&x,&b);CHKERRQ(ierr);
302   ierr = VecScatterCreate(b,shell->zrows,shell->zvals_w,shell->zrows,&shell->zvals_sct_r);CHKERRQ(ierr);
303 
304   /* col scatter: from right vector to left vector */
305   ierr = ISGetIndices(shell->zrows,&ridxs);CHKERRQ(ierr);
306   ierr = ISGetLocalSize(shell->zrows,&nr);CHKERRQ(ierr);
307   ierr = PetscMalloc1(nr,&gidxs);CHKERRQ(ierr);
308   for (i = 0, cum  = 0; i < nr; i++) {
309     if (ridxs[i] >= mat->cmap->N) continue;
310     gidxs[cum] = ridxs[i];
311     cum++;
312   }
313   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)mat),cum,gidxs,PETSC_OWN_POINTER,&is1);CHKERRQ(ierr);
314   ierr = VecScatterCreate(x,is1,shell->zvals_w,is1,&shell->zvals_sct_c);CHKERRQ(ierr);
315   ierr = ISDestroy(&is1);CHKERRQ(ierr);
316   ierr = VecDestroy(&x);CHKERRQ(ierr);
317   ierr = VecDestroy(&b);CHKERRQ(ierr);
318 
319   /* Expand/create index set of zeroed columns */
320   if (rc) {
321     ierr = PetscMalloc1(nc,&idxs);CHKERRQ(ierr);
322     for (i = 0; i < nc; i++) idxs[i] = cols[i] + cst;
323     ierr = ISCreateGeneral(PETSC_COMM_SELF,nc,idxs,PETSC_OWN_POINTER,&is1);CHKERRQ(ierr);
324     ierr = ISSort(is1);CHKERRQ(ierr);
325     if (shell->zcols) {
326       ierr = ISSum(shell->zcols,is1,&is2);CHKERRQ(ierr);
327       ierr = ISDestroy(&shell->zcols);CHKERRQ(ierr);
328       ierr = ISDestroy(&is1);CHKERRQ(ierr);
329       shell->zcols = is2;
330     } else shell->zcols = is1;
331   }
332   PetscFunctionReturn(0);
333 }
334 
335 static PetscErrorCode MatZeroRows_Shell(Mat mat,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
336 {
337   PetscInt       nr, *lrows;
338   PetscErrorCode ierr;
339 
340   PetscFunctionBegin;
341   if (x && b) {
342     Vec          xt;
343     PetscScalar *vals;
344     PetscInt    *gcols,i,st,nl,nc;
345 
346     ierr = PetscMalloc1(n,&gcols);CHKERRQ(ierr);
347     for (i = 0, nc = 0; i < n; i++) if (rows[i] < mat->cmap->N) gcols[nc++] = rows[i];
348 
349     ierr = MatCreateVecs(mat,&xt,NULL);CHKERRQ(ierr);
350     ierr = VecCopy(x,xt);CHKERRQ(ierr);
351     ierr = PetscCalloc1(nc,&vals);CHKERRQ(ierr);
352     ierr = VecSetValues(xt,nc,gcols,vals,INSERT_VALUES);CHKERRQ(ierr); /* xt = [x1, 0] */
353     ierr = PetscFree(vals);CHKERRQ(ierr);
354     ierr = VecAssemblyBegin(xt);CHKERRQ(ierr);
355     ierr = VecAssemblyEnd(xt);CHKERRQ(ierr);
356     ierr = VecAYPX(xt,-1.0,x);CHKERRQ(ierr);                           /* xt = [0, x2] */
357 
358     ierr = VecGetOwnershipRange(xt,&st,NULL);CHKERRQ(ierr);
359     ierr = VecGetLocalSize(xt,&nl);CHKERRQ(ierr);
360     ierr = VecGetArray(xt,&vals);CHKERRQ(ierr);
361     for (i = 0; i < nl; i++) {
362       PetscInt g = i + st;
363       if (g > mat->rmap->N) continue;
364       if (PetscAbsScalar(vals[i]) == 0.0) continue;
365       ierr = VecSetValue(b,g,diag*vals[i],INSERT_VALUES);CHKERRQ(ierr);
366     }
367     ierr = VecRestoreArray(xt,&vals);CHKERRQ(ierr);
368     ierr = VecAssemblyBegin(b);CHKERRQ(ierr);
369     ierr = VecAssemblyEnd(b);CHKERRQ(ierr);                            /* b  = [b1, x2 * diag] */
370     ierr = VecDestroy(&xt);CHKERRQ(ierr);
371     ierr = PetscFree(gcols);CHKERRQ(ierr);
372   }
373   ierr = PetscLayoutMapLocal(mat->rmap,n,rows,&nr,&lrows,NULL);CHKERRQ(ierr);
374   ierr = MatZeroRowsColumns_Local_Shell(mat,nr,lrows,0,NULL,diag,PETSC_FALSE);CHKERRQ(ierr);
375   ierr = PetscFree(lrows);CHKERRQ(ierr);
376   PetscFunctionReturn(0);
377 }
378 
379 static PetscErrorCode MatZeroRowsColumns_Shell(Mat mat,PetscInt n,const PetscInt rowscols[],PetscScalar diag,Vec x,Vec b)
380 {
381   PetscInt       *lrows, *lcols;
382   PetscInt       nr, nc;
383   PetscBool      congruent;
384   PetscErrorCode ierr;
385 
386   PetscFunctionBegin;
387   if (x && b) {
388     Vec          xt, bt;
389     PetscScalar *vals;
390     PetscInt    *grows,*gcols,i,st,nl;
391 
392     ierr = PetscMalloc2(n,&grows,n,&gcols);CHKERRQ(ierr);
393     for (i = 0, nr = 0; i < n; i++) if (rowscols[i] < mat->rmap->N) grows[nr++] = rowscols[i];
394     for (i = 0, nc = 0; i < n; i++) if (rowscols[i] < mat->cmap->N) gcols[nc++] = rowscols[i];
395     ierr = PetscCalloc1(n,&vals);CHKERRQ(ierr);
396 
397     ierr = MatCreateVecs(mat,&xt,&bt);CHKERRQ(ierr);
398     ierr = VecCopy(x,xt);CHKERRQ(ierr);
399     ierr = VecSetValues(xt,nc,gcols,vals,INSERT_VALUES);CHKERRQ(ierr); /* xt = [x1, 0] */
400     ierr = VecAssemblyBegin(xt);CHKERRQ(ierr);
401     ierr = VecAssemblyEnd(xt);CHKERRQ(ierr);
402     ierr = VecAXPY(xt,-1.0,x);CHKERRQ(ierr);                           /* xt = [0, -x2] */
403     ierr = MatMult(mat,xt,bt);CHKERRQ(ierr);                           /* bt = [-A12*x2,-A22*x2] */
404     ierr = VecSetValues(bt,nr,grows,vals,INSERT_VALUES);CHKERRQ(ierr); /* bt = [-A12*x2,0] */
405     ierr = VecAssemblyBegin(bt);CHKERRQ(ierr);
406     ierr = VecAssemblyEnd(bt);CHKERRQ(ierr);
407     ierr = VecAXPY(b,1.0,bt);CHKERRQ(ierr);                            /* b  = [b1 - A12*x2, b2] */
408     ierr = VecSetValues(bt,nr,grows,vals,INSERT_VALUES);CHKERRQ(ierr); /* b  = [b1 - A12*x2, 0] */
409     ierr = VecAssemblyBegin(bt);CHKERRQ(ierr);
410     ierr = VecAssemblyEnd(bt);CHKERRQ(ierr);
411     ierr = PetscFree(vals);CHKERRQ(ierr);
412 
413     ierr = VecGetOwnershipRange(xt,&st,NULL);CHKERRQ(ierr);
414     ierr = VecGetLocalSize(xt,&nl);CHKERRQ(ierr);
415     ierr = VecGetArray(xt,&vals);CHKERRQ(ierr);
416     for (i = 0; i < nl; i++) {
417       PetscInt g = i + st;
418       if (g > mat->rmap->N) continue;
419       if (PetscAbsScalar(vals[i]) == 0.0) continue;
420       ierr = VecSetValue(b,g,-diag*vals[i],INSERT_VALUES);CHKERRQ(ierr);
421     }
422     ierr = VecRestoreArray(xt,&vals);CHKERRQ(ierr);
423     ierr = VecAssemblyBegin(b);CHKERRQ(ierr);
424     ierr = VecAssemblyEnd(b);CHKERRQ(ierr);                            /* b  = [b1 - A12*x2, x2 * diag] */
425     ierr = VecDestroy(&xt);CHKERRQ(ierr);
426     ierr = VecDestroy(&bt);CHKERRQ(ierr);
427     ierr = PetscFree2(grows,gcols);CHKERRQ(ierr);
428   }
429   ierr = PetscLayoutMapLocal(mat->rmap,n,rowscols,&nr,&lrows,NULL);CHKERRQ(ierr);
430   ierr = MatHasCongruentLayouts(mat,&congruent);CHKERRQ(ierr);
431   if (congruent) {
432     nc    = nr;
433     lcols = lrows;
434   } else { /* MatZeroRowsColumns implicitly assumes the rowscols indices are for a square matrix, here we handle a more general case */
435     PetscInt i,nt,*t;
436 
437     ierr = PetscMalloc1(n,&t);CHKERRQ(ierr);
438     for (i = 0, nt = 0; i < n; i++) if (rowscols[i] < mat->cmap->N) t[nt++] = rowscols[i];
439     ierr = PetscLayoutMapLocal(mat->cmap,nt,t,&nc,&lcols,NULL);CHKERRQ(ierr);
440     ierr = PetscFree(t);CHKERRQ(ierr);
441   }
442   ierr = MatZeroRowsColumns_Local_Shell(mat,nr,lrows,nc,lcols,diag,PETSC_TRUE);CHKERRQ(ierr);
443   if (!congruent) {
444     ierr = PetscFree(lcols);CHKERRQ(ierr);
445   }
446   ierr = PetscFree(lrows);CHKERRQ(ierr);
447   PetscFunctionReturn(0);
448 }
449 
450 PetscErrorCode MatDestroy_Shell(Mat mat)
451 {
452   PetscErrorCode ierr;
453   Mat_Shell      *shell = (Mat_Shell*)mat->data;
454 
455   PetscFunctionBegin;
456   if (shell->ops->destroy) {
457     ierr = (*shell->ops->destroy)(mat);CHKERRQ(ierr);
458   }
459   ierr = VecDestroy(&shell->left);CHKERRQ(ierr);
460   ierr = VecDestroy(&shell->right);CHKERRQ(ierr);
461   ierr = VecDestroy(&shell->dshift);CHKERRQ(ierr);
462   ierr = VecDestroy(&shell->left_work);CHKERRQ(ierr);
463   ierr = VecDestroy(&shell->right_work);CHKERRQ(ierr);
464   ierr = VecDestroy(&shell->left_add_work);CHKERRQ(ierr);
465   ierr = VecDestroy(&shell->right_add_work);CHKERRQ(ierr);
466   ierr = MatDestroy(&shell->axpy);CHKERRQ(ierr);
467   ierr = VecDestroy(&shell->zvals_w);CHKERRQ(ierr);
468   ierr = VecDestroy(&shell->zvals);CHKERRQ(ierr);
469   ierr = VecScatterDestroy(&shell->zvals_sct_c);CHKERRQ(ierr);
470   ierr = VecScatterDestroy(&shell->zvals_sct_r);CHKERRQ(ierr);
471   ierr = ISDestroy(&shell->zrows);CHKERRQ(ierr);
472   ierr = ISDestroy(&shell->zcols);CHKERRQ(ierr);
473   ierr = PetscFree(mat->data);CHKERRQ(ierr);
474   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatShellGetContext_C",NULL);CHKERRQ(ierr);
475   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatShellSetContext_C",NULL);CHKERRQ(ierr);
476   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatShellSetManageScalingShifts_C",NULL);CHKERRQ(ierr);
477   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatShellSetOperation_C",NULL);CHKERRQ(ierr);
478   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatShellGetOperation_C",NULL);CHKERRQ(ierr);
479   PetscFunctionReturn(0);
480 }
481 
482 PetscErrorCode MatCopy_Shell(Mat A,Mat B,MatStructure str)
483 {
484   Mat_Shell       *shellA = (Mat_Shell*)A->data,*shellB = (Mat_Shell*)B->data;
485   PetscErrorCode  ierr;
486   PetscBool       matflg;
487 
488   PetscFunctionBegin;
489   ierr = PetscObjectTypeCompare((PetscObject)B,MATSHELL,&matflg);CHKERRQ(ierr);
490   if (!matflg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_NOTSAMETYPE,"Output matrix must be a MATSHELL");
491 
492   ierr = PetscMemcpy(B->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
493   ierr = PetscMemcpy(shellB->ops,shellA->ops,sizeof(struct _MatShellOps));CHKERRQ(ierr);
494 
495   if (shellA->ops->copy) {
496     ierr = (*shellA->ops->copy)(A,B,str);CHKERRQ(ierr);
497   }
498   shellB->vscale = shellA->vscale;
499   shellB->vshift = shellA->vshift;
500   if (shellA->dshift) {
501     if (!shellB->dshift) {
502       ierr = VecDuplicate(shellA->dshift,&shellB->dshift);CHKERRQ(ierr);
503     }
504     ierr = VecCopy(shellA->dshift,shellB->dshift);CHKERRQ(ierr);
505   } else {
506     ierr = VecDestroy(&shellB->dshift);CHKERRQ(ierr);
507   }
508   if (shellA->left) {
509     if (!shellB->left) {
510       ierr = VecDuplicate(shellA->left,&shellB->left);CHKERRQ(ierr);
511     }
512     ierr = VecCopy(shellA->left,shellB->left);CHKERRQ(ierr);
513   } else {
514     ierr = VecDestroy(&shellB->left);CHKERRQ(ierr);
515   }
516   if (shellA->right) {
517     if (!shellB->right) {
518       ierr = VecDuplicate(shellA->right,&shellB->right);CHKERRQ(ierr);
519     }
520     ierr = VecCopy(shellA->right,shellB->right);CHKERRQ(ierr);
521   } else {
522     ierr = VecDestroy(&shellB->right);CHKERRQ(ierr);
523   }
524   ierr = MatDestroy(&shellB->axpy);CHKERRQ(ierr);
525   if (shellA->axpy) {
526     ierr                = PetscObjectReference((PetscObject)shellA->axpy);CHKERRQ(ierr);
527     shellB->axpy        = shellA->axpy;
528     shellB->axpy_vscale = shellA->axpy_vscale;
529   }
530   if (shellA->zrows) {
531     ierr = ISDuplicate(shellA->zrows,&shellB->zrows);CHKERRQ(ierr);
532     if (shellA->zcols) {
533       ierr = ISDuplicate(shellA->zcols,&shellB->zcols);CHKERRQ(ierr);
534     }
535     ierr = VecDuplicate(shellA->zvals,&shellB->zvals);CHKERRQ(ierr);
536     ierr = VecCopy(shellA->zvals,shellB->zvals);CHKERRQ(ierr);
537     ierr = VecDuplicate(shellA->zvals_w,&shellB->zvals_w);CHKERRQ(ierr);
538     ierr = PetscObjectReference((PetscObject)shellA->zvals_sct_r);CHKERRQ(ierr);
539     ierr = PetscObjectReference((PetscObject)shellA->zvals_sct_c);CHKERRQ(ierr);
540     shellB->zvals_sct_r = shellA->zvals_sct_r;
541     shellB->zvals_sct_c = shellA->zvals_sct_c;
542   }
543   PetscFunctionReturn(0);
544 }
545 
546 PetscErrorCode MatDuplicate_Shell(Mat mat,MatDuplicateOption op,Mat *M)
547 {
548   PetscErrorCode ierr;
549   void           *ctx;
550 
551   PetscFunctionBegin;
552   ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr);
553   ierr = MatCreateShell(PetscObjectComm((PetscObject)mat),mat->rmap->n,mat->cmap->n,mat->rmap->N,mat->cmap->N,ctx,M);CHKERRQ(ierr);
554   if (op != MAT_DO_NOT_COPY_VALUES) {
555     ierr = MatCopy(mat,*M,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
556   }
557   PetscFunctionReturn(0);
558 }
559 
560 PetscErrorCode MatMult_Shell(Mat A,Vec x,Vec y)
561 {
562   Mat_Shell        *shell = (Mat_Shell*)A->data;
563   PetscErrorCode   ierr;
564   Vec              xx;
565   PetscObjectState instate,outstate;
566 
567   PetscFunctionBegin;
568   if (!shell->ops->mult) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Have not provided a MatMult() for this MATSHELL");
569   ierr = MatShellPreZeroRight(A,x,&xx);CHKERRQ(ierr);
570   ierr = MatShellPreScaleRight(A,xx,&xx);CHKERRQ(ierr);
571   ierr = PetscObjectStateGet((PetscObject)y, &instate);CHKERRQ(ierr);
572   ierr = (*shell->ops->mult)(A,xx,y);CHKERRQ(ierr);
573   ierr = PetscObjectStateGet((PetscObject)y, &outstate);CHKERRQ(ierr);
574   if (instate == outstate) {
575     /* increase the state of the output vector since the user did not update its state themself as should have been done */
576     ierr = PetscObjectStateIncrease((PetscObject)y);CHKERRQ(ierr);
577   }
578   ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr);
579   ierr = MatShellPostScaleLeft(A,y);CHKERRQ(ierr);
580   ierr = MatShellPostZeroLeft(A,y);CHKERRQ(ierr);
581 
582   if (shell->axpy) {
583     if (!shell->left_work) {ierr = MatCreateVecs(A,&shell->left_work,NULL);CHKERRQ(ierr);}
584     ierr = MatMult(shell->axpy,x,shell->left_work);CHKERRQ(ierr);
585     ierr = VecAXPY(y,shell->axpy_vscale,shell->left_work);CHKERRQ(ierr);
586   }
587   PetscFunctionReturn(0);
588 }
589 
590 PetscErrorCode MatMultAdd_Shell(Mat A,Vec x,Vec y,Vec z)
591 {
592   Mat_Shell      *shell = (Mat_Shell*)A->data;
593   PetscErrorCode ierr;
594 
595   PetscFunctionBegin;
596   if (y == z) {
597     if (!shell->right_add_work) {ierr = VecDuplicate(z,&shell->right_add_work);CHKERRQ(ierr);}
598     ierr = MatMult(A,x,shell->right_add_work);CHKERRQ(ierr);
599     ierr = VecAXPY(z,1.0,shell->right_add_work);CHKERRQ(ierr);
600   } else {
601     ierr = MatMult(A,x,z);CHKERRQ(ierr);
602     ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr);
603   }
604   PetscFunctionReturn(0);
605 }
606 
607 PetscErrorCode MatMultTranspose_Shell(Mat A,Vec x,Vec y)
608 {
609   Mat_Shell        *shell = (Mat_Shell*)A->data;
610   PetscErrorCode   ierr;
611   Vec              xx;
612   PetscObjectState instate,outstate;
613 
614   PetscFunctionBegin;
615   ierr = MatShellPreZeroLeft(A,x,&xx);CHKERRQ(ierr);
616   ierr = MatShellPreScaleLeft(A,xx,&xx);CHKERRQ(ierr);
617   ierr = PetscObjectStateGet((PetscObject)y, &instate);CHKERRQ(ierr);
618   if (!shell->ops->multtranspose) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Have not provided a MatMultTranspose() for this MATSHELL");
619   ierr = (*shell->ops->multtranspose)(A,xx,y);CHKERRQ(ierr);
620   ierr = PetscObjectStateGet((PetscObject)y, &outstate);CHKERRQ(ierr);
621   if (instate == outstate) {
622     /* increase the state of the output vector since the user did not update its state themself as should have been done */
623     ierr = PetscObjectStateIncrease((PetscObject)y);CHKERRQ(ierr);
624   }
625   ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr);
626   ierr = MatShellPostScaleRight(A,y);CHKERRQ(ierr);
627   ierr = MatShellPostZeroRight(A,y);CHKERRQ(ierr);
628 
629   if (shell->axpy) {
630     if (!shell->right_work) {ierr = MatCreateVecs(A,NULL,&shell->right_work);CHKERRQ(ierr);}
631     ierr = MatMultTranspose(shell->axpy,x,shell->right_work);CHKERRQ(ierr);
632     ierr = VecAXPY(y,shell->axpy_vscale,shell->right_work);CHKERRQ(ierr);
633   }
634   PetscFunctionReturn(0);
635 }
636 
637 PetscErrorCode MatMultTransposeAdd_Shell(Mat A,Vec x,Vec y,Vec z)
638 {
639   Mat_Shell      *shell = (Mat_Shell*)A->data;
640   PetscErrorCode ierr;
641 
642   PetscFunctionBegin;
643   if (y == z) {
644     if (!shell->left_add_work) {ierr = VecDuplicate(z,&shell->left_add_work);CHKERRQ(ierr);}
645     ierr = MatMultTranspose(A,x,shell->left_add_work);CHKERRQ(ierr);
646     ierr = VecAXPY(z,1.0,shell->left_add_work);CHKERRQ(ierr);
647   } else {
648     ierr = MatMultTranspose(A,x,z);CHKERRQ(ierr);
649     ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr);
650   }
651   PetscFunctionReturn(0);
652 }
653 
654 /*
655           diag(left)(vscale*A + diag(dshift) + vshift I)diag(right)
656 */
657 PetscErrorCode MatGetDiagonal_Shell(Mat A,Vec v)
658 {
659   Mat_Shell      *shell = (Mat_Shell*)A->data;
660   PetscErrorCode ierr;
661 
662   PetscFunctionBegin;
663   if (shell->ops->getdiagonal) {
664     ierr = (*shell->ops->getdiagonal)(A,v);CHKERRQ(ierr);
665   } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Must provide shell matrix with routine to return diagonal using\nMatShellSetOperation(S,MATOP_GET_DIAGONAL,...)");
666   ierr = VecScale(v,shell->vscale);CHKERRQ(ierr);
667   if (shell->dshift) {
668     ierr = VecAXPY(v,1.0,shell->dshift);CHKERRQ(ierr);
669   }
670   ierr = VecShift(v,shell->vshift);CHKERRQ(ierr);
671   if (shell->left)  {ierr = VecPointwiseMult(v,v,shell->left);CHKERRQ(ierr);}
672   if (shell->right) {ierr = VecPointwiseMult(v,v,shell->right);CHKERRQ(ierr);}
673   if (shell->zrows) {
674     ierr = VecScatterBegin(shell->zvals_sct_r,shell->zvals,v,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
675     ierr = VecScatterEnd(shell->zvals_sct_r,shell->zvals,v,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
676   }
677   if (shell->axpy) {
678     if (!shell->left_work) {ierr = VecDuplicate(v,&shell->left_work);CHKERRQ(ierr);}
679     ierr = MatGetDiagonal(shell->axpy,shell->left_work);CHKERRQ(ierr);
680     ierr = VecAXPY(v,shell->axpy_vscale,shell->left_work);CHKERRQ(ierr);
681   }
682   PetscFunctionReturn(0);
683 }
684 
685 PetscErrorCode MatShift_Shell(Mat Y,PetscScalar a)
686 {
687   Mat_Shell      *shell = (Mat_Shell*)Y->data;
688   PetscErrorCode ierr;
689   PetscBool      flg;
690 
691   PetscFunctionBegin;
692   ierr = MatHasCongruentLayouts(Y,&flg);CHKERRQ(ierr);
693   if (!flg) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_SUP,"Cannot shift shell matrix if it is not congruent");
694   if (shell->left || shell->right) {
695     if (!shell->dshift) {
696       ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->dshift);CHKERRQ(ierr);
697       ierr = VecSet(shell->dshift,a);CHKERRQ(ierr);
698     } else {
699       if (shell->left)  {ierr = VecPointwiseMult(shell->dshift,shell->dshift,shell->left);CHKERRQ(ierr);}
700       if (shell->right) {ierr = VecPointwiseMult(shell->dshift,shell->dshift,shell->right);CHKERRQ(ierr);}
701       ierr = VecShift(shell->dshift,a);CHKERRQ(ierr);
702     }
703     if (shell->left)  {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->left);CHKERRQ(ierr);}
704     if (shell->right) {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->right);CHKERRQ(ierr);}
705   } else shell->vshift += a;
706   if (shell->zrows) {
707     ierr = VecShift(shell->zvals,a);CHKERRQ(ierr);
708   }
709   PetscFunctionReturn(0);
710 }
711 
712 PetscErrorCode MatDiagonalSet_Shell_Private(Mat A,Vec D,PetscScalar s)
713 {
714   Mat_Shell      *shell = (Mat_Shell*)A->data;
715   PetscErrorCode ierr;
716 
717   PetscFunctionBegin;
718   if (!shell->dshift) {ierr = VecDuplicate(D,&shell->dshift);CHKERRQ(ierr);}
719   if (shell->left || shell->right) {
720     if (!shell->right_work) {ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->right_work);CHKERRQ(ierr);}
721     if (shell->left && shell->right)  {
722       ierr = VecPointwiseDivide(shell->right_work,D,shell->left);CHKERRQ(ierr);
723       ierr = VecPointwiseDivide(shell->right_work,shell->right_work,shell->right);CHKERRQ(ierr);
724     } else if (shell->left) {
725       ierr = VecPointwiseDivide(shell->right_work,D,shell->left);CHKERRQ(ierr);
726     } else {
727       ierr = VecPointwiseDivide(shell->right_work,D,shell->right);CHKERRQ(ierr);
728     }
729     ierr = VecAXPY(shell->dshift,s,shell->right_work);CHKERRQ(ierr);
730   } else {
731     ierr = VecAXPY(shell->dshift,s,D);CHKERRQ(ierr);
732   }
733   PetscFunctionReturn(0);
734 }
735 
736 PetscErrorCode MatDiagonalSet_Shell(Mat A,Vec D,InsertMode ins)
737 {
738   Mat_Shell      *shell = (Mat_Shell*)A->data;
739   Vec            d;
740   PetscErrorCode ierr;
741   PetscBool      flg;
742 
743   PetscFunctionBegin;
744   ierr = MatHasCongruentLayouts(A,&flg);CHKERRQ(ierr);
745   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot diagonal set or shift shell matrix if it is not congruent");
746   if (ins == INSERT_VALUES) {
747     if (!A->ops->getdiagonal) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Operation MATOP_GETDIAGONAL must be set first");
748     ierr = VecDuplicate(D,&d);CHKERRQ(ierr);
749     ierr = MatGetDiagonal(A,d);CHKERRQ(ierr);
750     ierr = MatDiagonalSet_Shell_Private(A,d,-1.);CHKERRQ(ierr);
751     ierr = MatDiagonalSet_Shell_Private(A,D,1.);CHKERRQ(ierr);
752     ierr = VecDestroy(&d);CHKERRQ(ierr);
753     if (shell->zrows) {
754       ierr = VecCopy(D,shell->zvals);CHKERRQ(ierr);
755     }
756   } else {
757     ierr = MatDiagonalSet_Shell_Private(A,D,1.);CHKERRQ(ierr);
758     if (shell->zrows) {
759       ierr = VecAXPY(shell->zvals,1.0,D);CHKERRQ(ierr);
760     }
761   }
762   PetscFunctionReturn(0);
763 }
764 
765 PetscErrorCode MatScale_Shell(Mat Y,PetscScalar a)
766 {
767   Mat_Shell      *shell = (Mat_Shell*)Y->data;
768   PetscErrorCode ierr;
769 
770   PetscFunctionBegin;
771   shell->vscale *= a;
772   shell->vshift *= a;
773   if (shell->dshift) {
774     ierr = VecScale(shell->dshift,a);CHKERRQ(ierr);
775   }
776   shell->axpy_vscale *= a;
777   if (shell->zrows) {
778     ierr = VecScale(shell->zvals,a);CHKERRQ(ierr);
779   }
780   PetscFunctionReturn(0);
781 }
782 
783 static PetscErrorCode MatDiagonalScale_Shell(Mat Y,Vec left,Vec right)
784 {
785   Mat_Shell      *shell = (Mat_Shell*)Y->data;
786   PetscErrorCode ierr;
787 
788   PetscFunctionBegin;
789   if (shell->axpy) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_SUP,"Cannot diagonal scale MATSHELL after MatAXPY operation");
790   if (left) {
791     if (!shell->left) {
792       ierr = VecDuplicate(left,&shell->left);CHKERRQ(ierr);
793       ierr = VecCopy(left,shell->left);CHKERRQ(ierr);
794     } else {
795       ierr = VecPointwiseMult(shell->left,shell->left,left);CHKERRQ(ierr);
796     }
797     if (shell->zrows) {
798       ierr = VecPointwiseMult(shell->zvals,shell->zvals,left);CHKERRQ(ierr);
799     }
800   }
801   if (right) {
802     if (!shell->right) {
803       ierr = VecDuplicate(right,&shell->right);CHKERRQ(ierr);
804       ierr = VecCopy(right,shell->right);CHKERRQ(ierr);
805     } else {
806       ierr = VecPointwiseMult(shell->right,shell->right,right);CHKERRQ(ierr);
807     }
808     if (shell->zrows) {
809       if (!shell->left_work) {
810         ierr = MatCreateVecs(Y,NULL,&shell->left_work);CHKERRQ(ierr);
811       }
812       ierr = VecSet(shell->zvals_w,1.0);CHKERRQ(ierr);
813       ierr = VecScatterBegin(shell->zvals_sct_c,right,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
814       ierr = VecScatterEnd(shell->zvals_sct_c,right,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
815       ierr = VecPointwiseMult(shell->zvals,shell->zvals,shell->zvals_w);CHKERRQ(ierr);
816     }
817   }
818   PetscFunctionReturn(0);
819 }
820 
821 PetscErrorCode MatAssemblyEnd_Shell(Mat Y,MatAssemblyType t)
822 {
823   Mat_Shell      *shell = (Mat_Shell*)Y->data;
824   PetscErrorCode ierr;
825 
826   PetscFunctionBegin;
827   if (t == MAT_FINAL_ASSEMBLY) {
828     shell->vshift = 0.0;
829     shell->vscale = 1.0;
830     ierr = VecDestroy(&shell->dshift);CHKERRQ(ierr);
831     ierr = VecDestroy(&shell->left);CHKERRQ(ierr);
832     ierr = VecDestroy(&shell->right);CHKERRQ(ierr);
833     ierr = MatDestroy(&shell->axpy);CHKERRQ(ierr);
834     ierr = VecScatterDestroy(&shell->zvals_sct_c);CHKERRQ(ierr);
835     ierr = VecScatterDestroy(&shell->zvals_sct_r);CHKERRQ(ierr);
836     ierr = ISDestroy(&shell->zrows);CHKERRQ(ierr);
837     ierr = ISDestroy(&shell->zcols);CHKERRQ(ierr);
838   }
839   PetscFunctionReturn(0);
840 }
841 
842 static PetscErrorCode MatMissingDiagonal_Shell(Mat A,PetscBool  *missing,PetscInt *d)
843 {
844   PetscFunctionBegin;
845   *missing = PETSC_FALSE;
846   PetscFunctionReturn(0);
847 }
848 
849 PetscErrorCode MatAXPY_Shell(Mat Y,PetscScalar a,Mat X,MatStructure str)
850 {
851   Mat_Shell      *shell = (Mat_Shell*)Y->data;
852   PetscErrorCode ierr;
853 
854   PetscFunctionBegin;
855   ierr = PetscObjectReference((PetscObject)X);CHKERRQ(ierr);
856   ierr = MatDestroy(&shell->axpy);CHKERRQ(ierr);
857   shell->axpy        = X;
858   shell->axpy_vscale = a;
859   PetscFunctionReturn(0);
860 }
861 
862 static struct _MatOps MatOps_Values = {0,
863                                        0,
864                                        0,
865                                        0,
866                                 /* 4*/ MatMultAdd_Shell,
867                                        0,
868                                        MatMultTransposeAdd_Shell,
869                                        0,
870                                        0,
871                                        0,
872                                 /*10*/ 0,
873                                        0,
874                                        0,
875                                        0,
876                                        0,
877                                 /*15*/ 0,
878                                        0,
879                                        0,
880                                        MatDiagonalScale_Shell,
881                                        0,
882                                 /*20*/ 0,
883                                        MatAssemblyEnd_Shell,
884                                        0,
885                                        0,
886                                 /*24*/ MatZeroRows_Shell,
887                                        0,
888                                        0,
889                                        0,
890                                        0,
891                                 /*29*/ 0,
892                                        0,
893                                        0,
894                                        0,
895                                        0,
896                                 /*34*/ MatDuplicate_Shell,
897                                        0,
898                                        0,
899                                        0,
900                                        0,
901                                 /*39*/ MatAXPY_Shell,
902                                        0,
903                                        0,
904                                        0,
905                                        MatCopy_Shell,
906                                 /*44*/ 0,
907                                        MatScale_Shell,
908                                        MatShift_Shell,
909                                        MatDiagonalSet_Shell,
910                                        MatZeroRowsColumns_Shell,
911                                 /*49*/ 0,
912                                        0,
913                                        0,
914                                        0,
915                                        0,
916                                 /*54*/ 0,
917                                        0,
918                                        0,
919                                        0,
920                                        0,
921                                 /*59*/ 0,
922                                        MatDestroy_Shell,
923                                        0,
924                                        MatConvertFrom_Shell,
925                                        0,
926                                 /*64*/ 0,
927                                        0,
928                                        0,
929                                        0,
930                                        0,
931                                 /*69*/ 0,
932                                        0,
933                                        MatConvert_Shell,
934                                        0,
935                                        0,
936                                 /*74*/ 0,
937                                        0,
938                                        0,
939                                        0,
940                                        0,
941                                 /*79*/ 0,
942                                        0,
943                                        0,
944                                        0,
945                                        0,
946                                 /*84*/ 0,
947                                        0,
948                                        0,
949                                        0,
950                                        0,
951                                 /*89*/ 0,
952                                        0,
953                                        0,
954                                        0,
955                                        0,
956                                 /*94*/ 0,
957                                        0,
958                                        0,
959                                        0,
960                                        0,
961                                 /*99*/ 0,
962                                        0,
963                                        0,
964                                        0,
965                                        0,
966                                /*104*/ 0,
967                                        0,
968                                        0,
969                                        0,
970                                        0,
971                                /*109*/ 0,
972                                        0,
973                                        0,
974                                        0,
975                                        MatMissingDiagonal_Shell,
976                                /*114*/ 0,
977                                        0,
978                                        0,
979                                        0,
980                                        0,
981                                /*119*/ 0,
982                                        0,
983                                        0,
984                                        0,
985                                        0,
986                                /*124*/ 0,
987                                        0,
988                                        0,
989                                        0,
990                                        0,
991                                /*129*/ 0,
992                                        0,
993                                        0,
994                                        0,
995                                        0,
996                                /*134*/ 0,
997                                        0,
998                                        0,
999                                        0,
1000                                        0,
1001                                /*139*/ 0,
1002                                        0,
1003                                        0
1004 };
1005 
1006 PetscErrorCode  MatShellSetContext_Shell(Mat mat,void *ctx)
1007 {
1008   Mat_Shell      *shell = (Mat_Shell*)mat->data;
1009 
1010   PetscFunctionBegin;
1011   shell->ctx = ctx;
1012   PetscFunctionReturn(0);
1013 }
1014 
1015 PetscErrorCode MatShellSetManageScalingShifts_Shell(Mat A)
1016 {
1017   Mat_Shell      *shell = (Mat_Shell*)A->data;
1018 
1019   PetscFunctionBegin;
1020   shell->managescalingshifts = PETSC_FALSE;
1021   A->ops->diagonalset   = NULL;
1022   A->ops->diagonalscale = NULL;
1023   A->ops->scale         = NULL;
1024   A->ops->shift         = NULL;
1025   A->ops->axpy          = NULL;
1026   PetscFunctionReturn(0);
1027 }
1028 
1029 PetscErrorCode MatShellSetOperation_Shell(Mat mat,MatOperation op,void (*f)(void))
1030 {
1031   Mat_Shell      *shell = (Mat_Shell*)mat->data;
1032 
1033   PetscFunctionBegin;
1034   switch (op) {
1035   case MATOP_DESTROY:
1036     shell->ops->destroy = (PetscErrorCode (*)(Mat))f;
1037     break;
1038   case MATOP_VIEW:
1039     if (!mat->ops->viewnative) {
1040       mat->ops->viewnative = mat->ops->view;
1041     }
1042     mat->ops->view = (PetscErrorCode (*)(Mat,PetscViewer))f;
1043     break;
1044   case MATOP_COPY:
1045     shell->ops->copy = (PetscErrorCode (*)(Mat,Mat,MatStructure))f;
1046     break;
1047   case MATOP_DIAGONAL_SET:
1048   case MATOP_DIAGONAL_SCALE:
1049   case MATOP_SHIFT:
1050   case MATOP_SCALE:
1051   case MATOP_AXPY:
1052   case MATOP_ZERO_ROWS:
1053   case MATOP_ZERO_ROWS_COLUMNS:
1054     if (shell->managescalingshifts) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"MATSHELL is managing scalings and shifts, see MatShellSetManageScalingShifts()");
1055     (((void(**)(void))mat->ops)[op]) = f;
1056     break;
1057   case MATOP_GET_DIAGONAL:
1058     if (shell->managescalingshifts) {
1059       shell->ops->getdiagonal = (PetscErrorCode (*)(Mat,Vec))f;
1060       mat->ops->getdiagonal   = MatGetDiagonal_Shell;
1061     } else {
1062       shell->ops->getdiagonal = NULL;
1063       mat->ops->getdiagonal   = (PetscErrorCode (*)(Mat,Vec))f;
1064     }
1065     break;
1066   case MATOP_MULT:
1067     if (shell->managescalingshifts) {
1068       shell->ops->mult = (PetscErrorCode (*)(Mat,Vec,Vec))f;
1069       mat->ops->mult   = MatMult_Shell;
1070     } else {
1071       shell->ops->mult = NULL;
1072       mat->ops->mult   = (PetscErrorCode (*)(Mat,Vec,Vec))f;
1073     }
1074     break;
1075   case MATOP_MULT_TRANSPOSE:
1076     if (shell->managescalingshifts) {
1077       shell->ops->multtranspose = (PetscErrorCode (*)(Mat,Vec,Vec))f;
1078       mat->ops->multtranspose   = MatMultTranspose_Shell;
1079     } else {
1080       shell->ops->multtranspose = NULL;
1081       mat->ops->multtranspose   = (PetscErrorCode (*)(Mat,Vec,Vec))f;
1082     }
1083     break;
1084   default:
1085     (((void(**)(void))mat->ops)[op]) = f;
1086     break;
1087   }
1088   PetscFunctionReturn(0);
1089 }
1090 
1091 PetscErrorCode MatShellGetOperation_Shell(Mat mat,MatOperation op,void(**f)(void))
1092 {
1093   Mat_Shell      *shell = (Mat_Shell*)mat->data;
1094 
1095   PetscFunctionBegin;
1096   switch (op) {
1097   case MATOP_DESTROY:
1098     *f = (void (*)(void))shell->ops->destroy;
1099     break;
1100   case MATOP_VIEW:
1101     *f = (void (*)(void))mat->ops->view;
1102     break;
1103   case MATOP_COPY:
1104     *f = (void (*)(void))shell->ops->copy;
1105     break;
1106   case MATOP_DIAGONAL_SET:
1107   case MATOP_DIAGONAL_SCALE:
1108   case MATOP_SHIFT:
1109   case MATOP_SCALE:
1110   case MATOP_AXPY:
1111   case MATOP_ZERO_ROWS:
1112   case MATOP_ZERO_ROWS_COLUMNS:
1113     *f = (((void (**)(void))mat->ops)[op]);
1114     break;
1115   case MATOP_GET_DIAGONAL:
1116     if (shell->ops->getdiagonal)
1117       *f = (void (*)(void))shell->ops->getdiagonal;
1118     else
1119       *f = (((void (**)(void))mat->ops)[op]);
1120     break;
1121   case MATOP_MULT:
1122     if (shell->ops->mult)
1123       *f = (void (*)(void))shell->ops->mult;
1124     else
1125       *f = (((void (**)(void))mat->ops)[op]);
1126     break;
1127   case MATOP_MULT_TRANSPOSE:
1128     if (shell->ops->multtranspose)
1129       *f = (void (*)(void))shell->ops->multtranspose;
1130     else
1131       *f = (((void (**)(void))mat->ops)[op]);
1132     break;
1133   default:
1134     *f = (((void (**)(void))mat->ops)[op]);
1135   }
1136   PetscFunctionReturn(0);
1137 }
1138 
1139 /*MC
1140    MATSHELL - MATSHELL = "shell" - A matrix type to be used to define your own matrix type -- perhaps matrix free.
1141 
1142   Level: advanced
1143 
1144 .seealso: MatCreateShell()
1145 M*/
1146 
1147 PETSC_EXTERN PetscErrorCode MatCreate_Shell(Mat A)
1148 {
1149   Mat_Shell      *b;
1150   PetscErrorCode ierr;
1151 
1152   PetscFunctionBegin;
1153   ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1154 
1155   ierr    = PetscNewLog(A,&b);CHKERRQ(ierr);
1156   A->data = (void*)b;
1157 
1158   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
1159   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
1160 
1161   b->ctx                 = 0;
1162   b->vshift              = 0.0;
1163   b->vscale              = 1.0;
1164   b->managescalingshifts = PETSC_TRUE;
1165   A->assembled           = PETSC_TRUE;
1166   A->preallocated        = PETSC_FALSE;
1167 
1168   ierr = PetscObjectComposeFunction((PetscObject)A,"MatShellGetContext_C",MatShellGetContext_Shell);CHKERRQ(ierr);
1169   ierr = PetscObjectComposeFunction((PetscObject)A,"MatShellSetContext_C",MatShellSetContext_Shell);CHKERRQ(ierr);
1170   ierr = PetscObjectComposeFunction((PetscObject)A,"MatShellSetManageScalingShifts_C",MatShellSetManageScalingShifts_Shell);CHKERRQ(ierr);
1171   ierr = PetscObjectComposeFunction((PetscObject)A,"MatShellSetOperation_C",MatShellSetOperation_Shell);CHKERRQ(ierr);
1172   ierr = PetscObjectComposeFunction((PetscObject)A,"MatShellGetOperation_C",MatShellGetOperation_Shell);CHKERRQ(ierr);
1173   ierr = PetscObjectChangeTypeName((PetscObject)A,MATSHELL);CHKERRQ(ierr);
1174   PetscFunctionReturn(0);
1175 }
1176 
1177 /*@C
1178    MatCreateShell - Creates a new matrix class for use with a user-defined
1179    private data storage format.
1180 
1181   Collective
1182 
1183    Input Parameters:
1184 +  comm - MPI communicator
1185 .  m - number of local rows (must be given)
1186 .  n - number of local columns (must be given)
1187 .  M - number of global rows (may be PETSC_DETERMINE)
1188 .  N - number of global columns (may be PETSC_DETERMINE)
1189 -  ctx - pointer to data needed by the shell matrix routines
1190 
1191    Output Parameter:
1192 .  A - the matrix
1193 
1194    Level: advanced
1195 
1196   Usage:
1197 $    extern PetscErrorCode mult(Mat,Vec,Vec);
1198 $    MatCreateShell(comm,m,n,M,N,ctx,&mat);
1199 $    MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult);
1200 $    [ Use matrix for operations that have been set ]
1201 $    MatDestroy(mat);
1202 
1203    Notes:
1204    The shell matrix type is intended to provide a simple class to use
1205    with KSP (such as, for use with matrix-free methods). You should not
1206    use the shell type if you plan to define a complete matrix class.
1207 
1208    Fortran Notes:
1209     To use this from Fortran with a ctx you must write an interface definition for this
1210     function and for MatShellGetContext() that tells Fortran the Fortran derived data type you are passing
1211     in as the ctx argument.
1212 
1213    PETSc requires that matrices and vectors being used for certain
1214    operations are partitioned accordingly.  For example, when
1215    creating a shell matrix, A, that supports parallel matrix-vector
1216    products using MatMult(A,x,y) the user should set the number
1217    of local matrix rows to be the number of local elements of the
1218    corresponding result vector, y. Note that this is information is
1219    required for use of the matrix interface routines, even though
1220    the shell matrix may not actually be physically partitioned.
1221    For example,
1222 
1223 $
1224 $     Vec x, y
1225 $     extern PetscErrorCode mult(Mat,Vec,Vec);
1226 $     Mat A
1227 $
1228 $     VecCreateMPI(comm,PETSC_DECIDE,M,&y);
1229 $     VecCreateMPI(comm,PETSC_DECIDE,N,&x);
1230 $     VecGetLocalSize(y,&m);
1231 $     VecGetLocalSize(x,&n);
1232 $     MatCreateShell(comm,m,n,M,N,ctx,&A);
1233 $     MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult);
1234 $     MatMult(A,x,y);
1235 $     MatDestroy(&A);
1236 $     VecDestroy(&y);
1237 $     VecDestroy(&x);
1238 $
1239 
1240 
1241    MATSHELL handles MatShift(), MatDiagonalSet(), MatDiagonalScale(), MatAXPY(), MatScale(), MatZeroRows() and MatZeroRowsColumns() internally, so these
1242    operations cannot be overwritten unless MatShellSetManageScalingShifts() is called.
1243 
1244 
1245     For rectangular matrices do all the scalings and shifts make sense?
1246 
1247     Developers Notes:
1248     Regarding shifting and scaling. The general form is
1249 
1250           diag(left)(vscale*A + diag(dshift) + vshift I)diag(right)
1251 
1252       The order you apply the operations is important. For example if you have a dshift then
1253       apply a MatScale(s) you get s*vscale*A + s*diag(shift). But if you first scale and then shift
1254       you get s*vscale*A + diag(shift)
1255 
1256           A is the user provided function.
1257 
1258    KSP/PC uses changes in the Mat's "state" to decide if preconditioners need to be rebuilt: PCSetUp() only calls the setup() for
1259    for the PC implementation if the Mat state has increased from the previous call. Thus to get changes in a MATSHELL to trigger
1260    an update in the preconditioner you must call MatAssemblyBegin()/MatAssemblyEnd() or PetscObjectStateIncrease((PetscObject)mat);
1261    each time the MATSHELL matrix has changed.
1262 
1263    Calling MatAssemblyBegin()/MatAssemblyEnd() on a MATSHELL removes any previously supplied shift and scales that were provided
1264    with MatDiagonalSet(), MatShift(), MatScale(), or MatDiagonalScale().
1265 
1266 .seealso: MatShellSetOperation(), MatHasOperation(), MatShellGetContext(), MatShellSetContext(), MATSHELL, MatShellSetManageScalingShifts()
1267 @*/
1268 PetscErrorCode  MatCreateShell(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx,Mat *A)
1269 {
1270   PetscErrorCode ierr;
1271 
1272   PetscFunctionBegin;
1273   ierr = MatCreate(comm,A);CHKERRQ(ierr);
1274   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
1275   ierr = MatSetType(*A,MATSHELL);CHKERRQ(ierr);
1276   ierr = MatShellSetContext(*A,ctx);CHKERRQ(ierr);
1277   ierr = MatSetUp(*A);CHKERRQ(ierr);
1278   PetscFunctionReturn(0);
1279 }
1280 
1281 
1282 /*@
1283     MatShellSetContext - sets the context for a shell matrix
1284 
1285    Logically Collective on Mat
1286 
1287     Input Parameters:
1288 +   mat - the shell matrix
1289 -   ctx - the context
1290 
1291    Level: advanced
1292 
1293    Fortran Notes:
1294     To use this from Fortran you must write a Fortran interface definition for this
1295     function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.
1296 
1297 .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation()
1298 @*/
1299 PetscErrorCode  MatShellSetContext(Mat mat,void *ctx)
1300 {
1301   PetscErrorCode ierr;
1302 
1303   PetscFunctionBegin;
1304   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
1305   ierr = PetscUseMethod(mat,"MatShellSetContext_C",(Mat,void*),(mat,ctx));CHKERRQ(ierr);
1306   PetscFunctionReturn(0);
1307 }
1308 
1309 /*@
1310     MatShellSetManageScalingShifts - Allows the user to control the scaling and shift operations of the MATSHELL. Must be called immediately
1311           after MatCreateShell()
1312 
1313    Logically Collective on Mat
1314 
1315     Input Parameter:
1316 .   mat - the shell matrix
1317 
1318   Level: advanced
1319 
1320 .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext(), MatShellSetOperation()
1321 @*/
1322 PetscErrorCode MatShellSetManageScalingShifts(Mat A)
1323 {
1324   PetscErrorCode ierr;
1325 
1326   PetscFunctionBegin;
1327   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1328   ierr = PetscUseMethod(A,"MatShellSetManageScalingShifts_C",(Mat),(A));CHKERRQ(ierr);
1329   PetscFunctionReturn(0);
1330 }
1331 
1332 /*@C
1333     MatShellTestMult - Compares the multiply routine provided to the MATSHELL with differencing on a given function.
1334 
1335    Logically Collective on Mat
1336 
1337     Input Parameters:
1338 +   mat - the shell matrix
1339 .   f - the function
1340 .   base - differences are computed around this vector, see MatMFFDSetBase(), for Jacobians this is the point at which the Jacobian is being evaluated
1341 -   ctx - an optional context for the function
1342 
1343    Output Parameter:
1344 .   flg - PETSC_TRUE if the multiply is likely correct
1345 
1346    Options Database:
1347 .   -mat_shell_test_mult_view - print if any differences are detected between the products and print the difference
1348 
1349    Level: advanced
1350 
1351    Fortran Notes:
1352     Not supported from Fortran
1353 
1354 .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellTestMultTranspose()
1355 @*/
1356 PetscErrorCode  MatShellTestMult(Mat mat,PetscErrorCode (*f)(void*,Vec,Vec),Vec base,void *ctx,PetscBool *flg)
1357 {
1358   PetscErrorCode ierr;
1359   PetscInt       m,n;
1360   Mat            mf,Dmf,Dmat,Ddiff;
1361   PetscReal      Diffnorm,Dmfnorm;
1362   PetscBool      v = PETSC_FALSE, flag = PETSC_TRUE;
1363 
1364   PetscFunctionBegin;
1365   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
1366   ierr = PetscOptionsHasName(NULL,((PetscObject)mat)->prefix,"-mat_shell_test_mult_view",&v);CHKERRQ(ierr);
1367   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
1368   ierr = MatCreateMFFD(PetscObjectComm((PetscObject)mat),m,n,PETSC_DECIDE,PETSC_DECIDE,&mf);CHKERRQ(ierr);
1369   ierr = MatMFFDSetFunction(mf,f,ctx);CHKERRQ(ierr);
1370   ierr = MatMFFDSetBase(mf,base,NULL);CHKERRQ(ierr);
1371 
1372   ierr = MatComputeOperator(mf,MATAIJ,&Dmf);CHKERRQ(ierr);
1373   ierr = MatComputeOperator(mat,MATAIJ,&Dmat);CHKERRQ(ierr);
1374 
1375   ierr = MatDuplicate(Dmat,MAT_COPY_VALUES,&Ddiff);CHKERRQ(ierr);
1376   ierr = MatAXPY(Ddiff,-1.0,Dmf,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
1377   ierr = MatNorm(Ddiff,NORM_FROBENIUS,&Diffnorm);CHKERRQ(ierr);
1378   ierr = MatNorm(Dmf,NORM_FROBENIUS,&Dmfnorm);CHKERRQ(ierr);
1379   if (Diffnorm/Dmfnorm > 10*PETSC_SQRT_MACHINE_EPSILON) {
1380     flag = PETSC_FALSE;
1381     if (v) {
1382       ierr = PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL and matrix free multiple appear to produce different results.\n Norm Ratio %g Difference results followed by finite difference one\n",(double)(Diffnorm/Dmfnorm));CHKERRQ(ierr);
1383       ierr = MatViewFromOptions(Ddiff,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr);
1384       ierr = MatViewFromOptions(Dmf,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr);
1385       ierr = MatViewFromOptions(Dmat,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr);
1386     }
1387   } else if (v) {
1388     ierr = PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL and matrix free multiple appear to produce the same results\n");CHKERRQ(ierr);
1389   }
1390   if (flg) *flg = flag;
1391   ierr = MatDestroy(&Ddiff);CHKERRQ(ierr);
1392   ierr = MatDestroy(&mf);CHKERRQ(ierr);
1393   ierr = MatDestroy(&Dmf);CHKERRQ(ierr);
1394   ierr = MatDestroy(&Dmat);CHKERRQ(ierr);
1395   PetscFunctionReturn(0);
1396 }
1397 
1398 /*@C
1399     MatShellTestMultTranpose - Compares the multiply transpose routine provided to the MATSHELL with differencing on a given function.
1400 
1401    Logically Collective on Mat
1402 
1403     Input Parameters:
1404 +   mat - the shell matrix
1405 .   f - the function
1406 .   base - differences are computed around this vector, see MatMFFDSetBase(), for Jacobians this is the point at which the Jacobian is being evaluated
1407 -   ctx - an optional context for the function
1408 
1409    Output Parameter:
1410 .   flg - PETSC_TRUE if the multiply is likely correct
1411 
1412    Options Database:
1413 .   -mat_shell_test_mult_view - print if any differences are detected between the products and print the difference
1414 
1415    Level: advanced
1416 
1417    Fortran Notes:
1418     Not supported from Fortran
1419 
1420 .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellTestMult()
1421 @*/
1422 PetscErrorCode  MatShellTestMultTranspose(Mat mat,PetscErrorCode (*f)(void*,Vec,Vec),Vec base,void *ctx,PetscBool *flg)
1423 {
1424   PetscErrorCode ierr;
1425   Vec            x,y,z;
1426   PetscInt       m,n,M,N;
1427   Mat            mf,Dmf,Dmat,Ddiff;
1428   PetscReal      Diffnorm,Dmfnorm;
1429   PetscBool      v = PETSC_FALSE, flag = PETSC_TRUE;
1430 
1431   PetscFunctionBegin;
1432   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
1433   ierr = PetscOptionsHasName(NULL,((PetscObject)mat)->prefix,"-mat_shell_test_mult_transpose_view",&v);CHKERRQ(ierr);
1434   ierr = MatCreateVecs(mat,&x,&y);CHKERRQ(ierr);
1435   ierr = VecDuplicate(y,&z);CHKERRQ(ierr);
1436   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
1437   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
1438   ierr = MatCreateMFFD(PetscObjectComm((PetscObject)mat),m,n,M,N,&mf);CHKERRQ(ierr);
1439   ierr = MatMFFDSetFunction(mf,f,ctx);CHKERRQ(ierr);
1440   ierr = MatMFFDSetBase(mf,base,NULL);CHKERRQ(ierr);
1441   ierr = MatComputeOperator(mf,MATAIJ,&Dmf);CHKERRQ(ierr);
1442   ierr = MatTranspose(Dmf,MAT_INPLACE_MATRIX,&Dmf);CHKERRQ(ierr);
1443   ierr = MatComputeOperatorTranspose(mat,MATAIJ,&Dmat);CHKERRQ(ierr);
1444 
1445   ierr = MatDuplicate(Dmat,MAT_COPY_VALUES,&Ddiff);CHKERRQ(ierr);
1446   ierr = MatAXPY(Ddiff,-1.0,Dmf,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
1447   ierr = MatNorm(Ddiff,NORM_FROBENIUS,&Diffnorm);CHKERRQ(ierr);
1448   ierr = MatNorm(Dmf,NORM_FROBENIUS,&Dmfnorm);CHKERRQ(ierr);
1449   if (Diffnorm/Dmfnorm > 10*PETSC_SQRT_MACHINE_EPSILON) {
1450     flag = PETSC_FALSE;
1451     if (v) {
1452       ierr = PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL and matrix free multiple appear to produce different results.\n Norm Ratio %g Difference results followed by finite difference one\n",(double)(Diffnorm/Dmfnorm));CHKERRQ(ierr);
1453       ierr = MatViewFromOptions(Ddiff,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr);
1454       ierr = MatViewFromOptions(Dmf,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr);
1455       ierr = MatViewFromOptions(Dmat,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr);
1456     }
1457   } else if (v) {
1458     ierr = PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL transpose and matrix free multiple appear to produce the same results\n");CHKERRQ(ierr);
1459   }
1460   if (flg) *flg = flag;
1461   ierr = MatDestroy(&mf);CHKERRQ(ierr);
1462   ierr = MatDestroy(&Dmat);CHKERRQ(ierr);
1463   ierr = MatDestroy(&Ddiff);CHKERRQ(ierr);
1464   ierr = MatDestroy(&Dmf);CHKERRQ(ierr);
1465   ierr = VecDestroy(&x);CHKERRQ(ierr);
1466   ierr = VecDestroy(&y);CHKERRQ(ierr);
1467   ierr = VecDestroy(&z);CHKERRQ(ierr);
1468   PetscFunctionReturn(0);
1469 }
1470 
1471 /*@C
1472     MatShellSetOperation - Allows user to set a matrix operation for a shell matrix.
1473 
1474    Logically Collective on Mat
1475 
1476     Input Parameters:
1477 +   mat - the shell matrix
1478 .   op - the name of the operation
1479 -   g - the function that provides the operation.
1480 
1481    Level: advanced
1482 
1483     Usage:
1484 $      extern PetscErrorCode usermult(Mat,Vec,Vec);
1485 $      ierr = MatCreateShell(comm,m,n,M,N,ctx,&A);
1486 $      ierr = MatShellSetOperation(A,MATOP_MULT,(void(*)(void))usermult);
1487 
1488     Notes:
1489     See the file include/petscmat.h for a complete list of matrix
1490     operations, which all have the form MATOP_<OPERATION>, where
1491     <OPERATION> is the name (in all capital letters) of the
1492     user interface routine (e.g., MatMult() -> MATOP_MULT).
1493 
1494     All user-provided functions (except for MATOP_DESTROY) should have the same calling
1495     sequence as the usual matrix interface routines, since they
1496     are intended to be accessed via the usual matrix interface
1497     routines, e.g.,
1498 $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)
1499 
1500     In particular each function MUST return an error code of 0 on success and
1501     nonzero on failure.
1502 
1503     Within each user-defined routine, the user should call
1504     MatShellGetContext() to obtain the user-defined context that was
1505     set by MatCreateShell().
1506 
1507     Fortran Notes:
1508     For MatCreateVecs() the user code should check if the input left or right matrix is -1 and in that case not
1509        generate a matrix. See src/mat/examples/tests/ex120f.F
1510 
1511     Use MatSetOperation() to set an operation for any matrix type
1512 
1513 .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext(), MatSetOperation(), MatShellSetManageScalingShifts()
1514 @*/
1515 PetscErrorCode MatShellSetOperation(Mat mat,MatOperation op,void (*g)(void))
1516 {
1517   PetscErrorCode ierr;
1518 
1519   PetscFunctionBegin;
1520   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
1521   ierr = PetscUseMethod(mat,"MatShellSetOperation_C",(Mat,MatOperation,void (*)(void)),(mat,op,g));CHKERRQ(ierr);
1522   PetscFunctionReturn(0);
1523 }
1524 
1525 /*@C
1526     MatShellGetOperation - Gets a matrix function for a shell matrix.
1527 
1528     Not Collective
1529 
1530     Input Parameters:
1531 +   mat - the shell matrix
1532 -   op - the name of the operation
1533 
1534     Output Parameter:
1535 .   g - the function that provides the operation.
1536 
1537     Level: advanced
1538 
1539     Notes:
1540     See the file include/petscmat.h for a complete list of matrix
1541     operations, which all have the form MATOP_<OPERATION>, where
1542     <OPERATION> is the name (in all capital letters) of the
1543     user interface routine (e.g., MatMult() -> MATOP_MULT).
1544 
1545     All user-provided functions have the same calling
1546     sequence as the usual matrix interface routines, since they
1547     are intended to be accessed via the usual matrix interface
1548     routines, e.g.,
1549 $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)
1550 
1551     Within each user-defined routine, the user should call
1552     MatShellGetContext() to obtain the user-defined context that was
1553     set by MatCreateShell().
1554 
1555 .seealso: MatCreateShell(), MatShellGetContext(), MatShellSetOperation(), MatShellSetContext()
1556 @*/
1557 PetscErrorCode MatShellGetOperation(Mat mat,MatOperation op,void(**g)(void))
1558 {
1559   PetscErrorCode ierr;
1560 
1561   PetscFunctionBegin;
1562   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
1563   ierr = PetscUseMethod(mat,"MatShellGetOperation_C",(Mat,MatOperation,void (**)(void)),(mat,op,g));CHKERRQ(ierr);
1564   PetscFunctionReturn(0);
1565 }
1566