xref: /petsc/src/mat/impls/shell/shell.c (revision 0cc1a30777ba4e09cf41d097dacadf23dd8db198)
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 
10 struct _MatShellOps {
11   /*   0 */
12   PetscErrorCode (*mult)(Mat,Vec,Vec);
13   /*   5 */
14   PetscErrorCode (*multtranspose)(Mat,Vec,Vec);
15   /*  10 */
16   /*  15 */
17   PetscErrorCode (*getdiagonal)(Mat,Vec);
18   /*  20 */
19   /*  24 */
20   /*  29 */
21   /*  34 */
22   /*  39 */
23   PetscErrorCode (*axpy)(Mat,PetscScalar,Mat,MatStructure);
24   PetscErrorCode (*copy)(Mat,Mat,MatStructure);
25   /*  44 */
26   PetscErrorCode (*diagonalset)(Mat,Vec,InsertMode);
27   /*  49 */
28   /*  54 */
29   /*  59 */
30   PetscErrorCode (*destroy)(Mat);
31   /*  64 */
32   /*  69 */
33   /*  74 */
34   /*  79 */
35   /*  84 */
36   /*  89 */
37   /*  94 */
38   /*  99 */
39   /* 104 */
40   /* 109 */
41   /* 114 */
42   /* 119 */
43   /* 124 */
44   /* 129 */
45   /* 134 */
46   /* 139 */
47   /* 144 */
48 };
49 
50 typedef struct {
51   struct _MatShellOps ops[1];
52 
53   PetscScalar vscale,vshift;
54   Vec         dshift;
55   Vec         left,right;
56   Vec         left_work,right_work;
57   Vec         left_add_work,right_add_work;
58   Mat         axpy;
59   PetscScalar axpy_vscale;
60   PetscBool   managescalingshifts;                   /* The user will manage the scaling and shifts for the MATSHELL, not the default */
61   void        *ctx;
62 } Mat_Shell;
63 
64 /*
65       xx = diag(left)*x
66 */
67 static PetscErrorCode MatShellPreScaleLeft(Mat A,Vec x,Vec *xx)
68 {
69   Mat_Shell      *shell = (Mat_Shell*)A->data;
70   PetscErrorCode ierr;
71 
72   PetscFunctionBegin;
73   *xx = NULL;
74   if (!shell->left) {
75     *xx = x;
76   } else {
77     if (!shell->left_work) {ierr = VecDuplicate(shell->left,&shell->left_work);CHKERRQ(ierr);}
78     ierr = VecPointwiseMult(shell->left_work,x,shell->left);CHKERRQ(ierr);
79     *xx  = shell->left_work;
80   }
81   PetscFunctionReturn(0);
82 }
83 
84 /*
85      xx = diag(right)*x
86 */
87 static PetscErrorCode MatShellPreScaleRight(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->right) {
95     *xx = x;
96   } else {
97     if (!shell->right_work) {ierr = VecDuplicate(shell->right,&shell->right_work);CHKERRQ(ierr);}
98     ierr = VecPointwiseMult(shell->right_work,x,shell->right);CHKERRQ(ierr);
99     *xx  = shell->right_work;
100   }
101   PetscFunctionReturn(0);
102 }
103 
104 /*
105     x = diag(left)*x
106 */
107 static PetscErrorCode MatShellPostScaleLeft(Mat A,Vec x)
108 {
109   Mat_Shell      *shell = (Mat_Shell*)A->data;
110   PetscErrorCode ierr;
111 
112   PetscFunctionBegin;
113   if (shell->left) {ierr = VecPointwiseMult(x,x,shell->left);CHKERRQ(ierr);}
114   PetscFunctionReturn(0);
115 }
116 
117 /*
118     x = diag(right)*x
119 */
120 static PetscErrorCode MatShellPostScaleRight(Mat A,Vec x)
121 {
122   Mat_Shell      *shell = (Mat_Shell*)A->data;
123   PetscErrorCode ierr;
124 
125   PetscFunctionBegin;
126   if (shell->right) {ierr = VecPointwiseMult(x,x,shell->right);CHKERRQ(ierr);}
127   PetscFunctionReturn(0);
128 }
129 
130 /*
131          Y = vscale*Y + diag(dshift)*X + vshift*X
132 
133          On input Y already contains A*x
134 */
135 static PetscErrorCode MatShellShiftAndScale(Mat A,Vec X,Vec Y)
136 {
137   Mat_Shell      *shell = (Mat_Shell*)A->data;
138   PetscErrorCode ierr;
139 
140   PetscFunctionBegin;
141   if (shell->dshift) {          /* get arrays because there is no VecPointwiseMultAdd() */
142     PetscInt          i,m;
143     const PetscScalar *x,*d;
144     PetscScalar       *y;
145     ierr = VecGetLocalSize(X,&m);CHKERRQ(ierr);
146     ierr = VecGetArrayRead(shell->dshift,&d);CHKERRQ(ierr);
147     ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
148     ierr = VecGetArray(Y,&y);CHKERRQ(ierr);
149     for (i=0; i<m; i++) y[i] = shell->vscale*y[i] + d[i]*x[i];
150     ierr = VecRestoreArrayRead(shell->dshift,&d);CHKERRQ(ierr);
151     ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
152     ierr = VecRestoreArray(Y,&y);CHKERRQ(ierr);
153   } else {
154     ierr = VecScale(Y,shell->vscale);CHKERRQ(ierr);
155   }
156   if (shell->vshift != 0.0) {ierr = VecAXPY(Y,shell->vshift,X);CHKERRQ(ierr);} /* if test is for non-square matrices */
157   PetscFunctionReturn(0);
158 }
159 
160 /*@
161     MatShellGetContext - Returns the user-provided context associated with a shell matrix.
162 
163     Not Collective
164 
165     Input Parameter:
166 .   mat - the matrix, should have been created with MatCreateShell()
167 
168     Output Parameter:
169 .   ctx - the user provided context
170 
171     Level: advanced
172 
173    Fortran Notes: To use this from Fortran you must write a Fortran interface definition for this
174     function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.
175 
176 .keywords: matrix, shell, get, context
177 
178 .seealso: MatCreateShell(), MatShellSetOperation(), MatShellSetContext()
179 @*/
180 PetscErrorCode  MatShellGetContext(Mat mat,void *ctx)
181 {
182   PetscErrorCode ierr;
183   PetscBool      flg;
184 
185   PetscFunctionBegin;
186   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
187   PetscValidPointer(ctx,2);
188   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
189   if (flg) *(void**)ctx = ((Mat_Shell*)(mat->data))->ctx;
190   else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot get context from non-shell matrix");
191   PetscFunctionReturn(0);
192 }
193 
194 PetscErrorCode MatDestroy_Shell(Mat mat)
195 {
196   PetscErrorCode ierr;
197   Mat_Shell      *shell = (Mat_Shell*)mat->data;
198 
199   PetscFunctionBegin;
200   if (shell->ops->destroy) {
201     ierr = (*shell->ops->destroy)(mat);CHKERRQ(ierr);
202   }
203   ierr = VecDestroy(&shell->left);CHKERRQ(ierr);
204   ierr = VecDestroy(&shell->right);CHKERRQ(ierr);
205   ierr = VecDestroy(&shell->dshift);CHKERRQ(ierr);
206   ierr = VecDestroy(&shell->left_work);CHKERRQ(ierr);
207   ierr = VecDestroy(&shell->right_work);CHKERRQ(ierr);
208   ierr = VecDestroy(&shell->left_add_work);CHKERRQ(ierr);
209   ierr = VecDestroy(&shell->right_add_work);CHKERRQ(ierr);
210   ierr = MatDestroy(&shell->axpy);CHKERRQ(ierr);
211   ierr = PetscFree(mat->data);CHKERRQ(ierr);
212   PetscFunctionReturn(0);
213 }
214 
215 PetscErrorCode MatCopy_Shell(Mat A,Mat B,MatStructure str)
216 {
217   Mat_Shell       *shellA = (Mat_Shell*)A->data,*shellB = (Mat_Shell*)B->data;
218   PetscErrorCode  ierr;
219   PetscBool       matflg;
220 
221   PetscFunctionBegin;
222   ierr = PetscObjectTypeCompare((PetscObject)B,MATSHELL,&matflg);CHKERRQ(ierr);
223   if (!matflg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_NOTSAMETYPE,"Output matrix must be a MATSHELL");
224 
225   ierr = PetscMemcpy(B->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
226   ierr = PetscMemcpy(shellB->ops,shellA->ops,sizeof(struct _MatShellOps));CHKERRQ(ierr);
227 
228   if (shellA->ops->copy) {
229     ierr = (*shellA->ops->copy)(A,B,str);CHKERRQ(ierr);
230   }
231   shellB->vscale = shellA->vscale;
232   shellB->vshift = shellA->vshift;
233   if (shellA->dshift) {
234     if (!shellB->dshift) {
235       ierr = VecDuplicate(shellA->dshift,&shellB->dshift);CHKERRQ(ierr);
236     }
237     ierr = VecCopy(shellA->dshift,shellB->dshift);CHKERRQ(ierr);
238   } else {
239     ierr = VecDestroy(&shellB->dshift);CHKERRQ(ierr);
240   }
241   if (shellA->left) {
242     if (!shellB->left) {
243       ierr = VecDuplicate(shellA->left,&shellB->left);CHKERRQ(ierr);
244     }
245     ierr = VecCopy(shellA->left,shellB->left);CHKERRQ(ierr);
246   } else {
247     ierr = VecDestroy(&shellB->left);CHKERRQ(ierr);
248   }
249   if (shellA->right) {
250     if (!shellB->right) {
251       ierr = VecDuplicate(shellA->right,&shellB->right);CHKERRQ(ierr);
252     }
253     ierr = VecCopy(shellA->right,shellB->right);CHKERRQ(ierr);
254   } else {
255     ierr = VecDestroy(&shellB->right);CHKERRQ(ierr);
256   }
257   ierr = MatDestroy(&shellB->axpy);CHKERRQ(ierr);
258   if (shellA->axpy) {
259     ierr                 = PetscObjectReference((PetscObject)shellA->axpy);CHKERRQ(ierr);
260     shellB->axpy        = shellA->axpy;
261     shellB->axpy_vscale = shellA->axpy_vscale;
262   }
263   PetscFunctionReturn(0);
264 }
265 
266 PetscErrorCode MatDuplicate_Shell(Mat mat,MatDuplicateOption op,Mat *M)
267 {
268   PetscErrorCode ierr;
269   void           *ctx;
270 
271   PetscFunctionBegin;
272   ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr);
273   ierr = MatCreateShell(PetscObjectComm((PetscObject)mat),mat->rmap->n,mat->cmap->n,mat->rmap->N,mat->cmap->N,ctx,M);CHKERRQ(ierr);
274   ierr = MatCopy(mat,*M,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
275   PetscFunctionReturn(0);
276 }
277 
278 PetscErrorCode MatMult_Shell(Mat A,Vec x,Vec y)
279 {
280   Mat_Shell        *shell = (Mat_Shell*)A->data;
281   PetscErrorCode   ierr;
282   Vec              xx;
283   PetscObjectState instate,outstate;
284 
285   PetscFunctionBegin;
286   ierr = MatShellPreScaleRight(A,x,&xx);CHKERRQ(ierr);
287   ierr = PetscObjectStateGet((PetscObject)y, &instate);CHKERRQ(ierr);
288   ierr = (*shell->ops->mult)(A,xx,y);CHKERRQ(ierr);
289   ierr = PetscObjectStateGet((PetscObject)y, &outstate);CHKERRQ(ierr);
290   if (instate == outstate) {
291     /* increase the state of the output vector since the user did not update its state themself as should have been done */
292     ierr = PetscObjectStateIncrease((PetscObject)y);CHKERRQ(ierr);
293   }
294   ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr);
295   ierr = MatShellPostScaleLeft(A,y);CHKERRQ(ierr);
296 
297   if (shell->axpy) {
298     if (!shell->left_work) {ierr = MatCreateVecs(A,&shell->left_work,NULL);CHKERRQ(ierr);}
299     ierr = MatMult(shell->axpy,x,shell->left_work);CHKERRQ(ierr);
300     ierr = VecAXPY(y,shell->axpy_vscale,shell->left_work);CHKERRQ(ierr);
301   }
302   PetscFunctionReturn(0);
303 }
304 
305 PetscErrorCode MatMultAdd_Shell(Mat A,Vec x,Vec y,Vec z)
306 {
307   Mat_Shell      *shell = (Mat_Shell*)A->data;
308   PetscErrorCode ierr;
309 
310   PetscFunctionBegin;
311   if (y == z) {
312     if (!shell->right_add_work) {ierr = VecDuplicate(z,&shell->right_add_work);CHKERRQ(ierr);}
313     ierr = MatMult(A,x,shell->right_add_work);CHKERRQ(ierr);
314     ierr = VecAXPY(z,1.0,shell->right_add_work);CHKERRQ(ierr);
315   } else {
316     ierr = MatMult(A,x,z);CHKERRQ(ierr);
317     ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr);
318   }
319   PetscFunctionReturn(0);
320 }
321 
322 PetscErrorCode MatMultTranspose_Shell(Mat A,Vec x,Vec y)
323 {
324   Mat_Shell        *shell = (Mat_Shell*)A->data;
325   PetscErrorCode   ierr;
326   Vec              xx;
327   PetscObjectState instate,outstate;
328 
329   PetscFunctionBegin;
330   ierr = MatShellPreScaleLeft(A,x,&xx);CHKERRQ(ierr);
331   ierr = PetscObjectStateGet((PetscObject)y, &instate);CHKERRQ(ierr);
332   if (!shell->ops->multtranspose) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Have not provided a MatMultTranspose() for this MATSHELL");
333   ierr = (*shell->ops->multtranspose)(A,xx,y);CHKERRQ(ierr);
334   ierr = PetscObjectStateGet((PetscObject)y, &outstate);CHKERRQ(ierr);
335   if (instate == outstate) {
336     /* increase the state of the output vector since the user did not update its state themself as should have been done */
337     ierr = PetscObjectStateIncrease((PetscObject)y);CHKERRQ(ierr);
338   }
339   ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr);
340   ierr = MatShellPostScaleRight(A,y);CHKERRQ(ierr);
341   PetscFunctionReturn(0);
342 }
343 
344 PetscErrorCode MatMultTransposeAdd_Shell(Mat A,Vec x,Vec y,Vec z)
345 {
346   Mat_Shell      *shell = (Mat_Shell*)A->data;
347   PetscErrorCode ierr;
348 
349   PetscFunctionBegin;
350   if (y == z) {
351     if (!shell->left_add_work) {ierr = VecDuplicate(z,&shell->left_add_work);CHKERRQ(ierr);}
352     ierr = MatMultTranspose(A,x,shell->left_add_work);CHKERRQ(ierr);
353     ierr = VecWAXPY(z,1.0,shell->left_add_work,y);CHKERRQ(ierr);
354   } else {
355     ierr = MatMultTranspose(A,x,z);CHKERRQ(ierr);
356     ierr = VecAXPY(z,1.0,y);CHKERRQ(ierr);
357   }
358   PetscFunctionReturn(0);
359 }
360 
361 /*
362           diag(left)(vscale*A + diag(dshift) + vshift I)diag(right)
363 */
364 PetscErrorCode MatGetDiagonal_Shell(Mat A,Vec v)
365 {
366   Mat_Shell      *shell = (Mat_Shell*)A->data;
367   PetscErrorCode ierr;
368 
369   PetscFunctionBegin;
370   if (shell->ops->getdiagonal) {
371     ierr = (*shell->ops->getdiagonal)(A,v);CHKERRQ(ierr);
372   } else {
373     ierr = VecSet(v,0.0);CHKERRQ(ierr);
374   }
375   ierr = VecScale(v,shell->vscale);CHKERRQ(ierr);
376   if (shell->dshift) {
377     ierr = VecAXPY(v,1.0,shell->dshift);CHKERRQ(ierr);
378   }
379   ierr = VecShift(v,shell->vshift);CHKERRQ(ierr);
380   if (shell->left)  {ierr = VecPointwiseMult(v,v,shell->left);CHKERRQ(ierr);}
381   if (shell->right) {ierr = VecPointwiseMult(v,v,shell->right);CHKERRQ(ierr);}
382   if (shell->axpy) {
383     if (!shell->left_work) {ierr = VecDuplicate(v,&shell->left_work);CHKERRQ(ierr);}
384     ierr = MatGetDiagonal(shell->axpy,shell->left_work);CHKERRQ(ierr);
385     ierr = VecAXPY(v,shell->axpy_vscale,shell->left_work);CHKERRQ(ierr);
386   }
387   PetscFunctionReturn(0);
388 }
389 
390 PetscErrorCode MatShift_Shell(Mat Y,PetscScalar a)
391 {
392   Mat_Shell      *shell = (Mat_Shell*)Y->data;
393   PetscErrorCode ierr;
394 
395   PetscFunctionBegin;
396   if (shell->left || shell->right) {
397     if (!shell->dshift) {
398       ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->dshift);CHKERRQ(ierr);
399       ierr = VecSet(shell->dshift,a);CHKERRQ(ierr);
400     } else {
401       if (shell->left)  {ierr = VecPointwiseMult(shell->dshift,shell->dshift,shell->left);CHKERRQ(ierr);}
402       if (shell->right) {ierr = VecPointwiseMult(shell->dshift,shell->dshift,shell->right);CHKERRQ(ierr);}
403       ierr = VecShift(shell->dshift,a);CHKERRQ(ierr);
404     }
405     if (shell->left)  {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->left);CHKERRQ(ierr);}
406     if (shell->right) {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->right);CHKERRQ(ierr);}
407   } else shell->vshift += a;
408   PetscFunctionReturn(0);
409 }
410 
411 PetscErrorCode MatDiagonalSet_Shell(Mat A,Vec D,InsertMode ins)
412 {
413   Mat_Shell      *shell = (Mat_Shell*)A->data;
414   PetscErrorCode ierr;
415 
416   PetscFunctionBegin;
417   if (ins == INSERT_VALUES) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE, "Operation not supported with INSERT_VALUES");
418   if (!shell->dshift) {ierr = VecDuplicate(D,&shell->dshift);CHKERRQ(ierr);}
419   if (shell->left || shell->right) {
420     if (!shell->right_work) {ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->right_work);CHKERRQ(ierr);}
421     if (shell->left && shell->right)  {
422       ierr = VecPointwiseDivide(shell->right_work,D,shell->left);CHKERRQ(ierr);
423       ierr = VecPointwiseDivide(shell->right_work,shell->right_work,shell->right);CHKERRQ(ierr);
424     } else if (shell->left) {
425       ierr = VecPointwiseDivide(shell->right_work,D,shell->left);CHKERRQ(ierr);
426     } else {
427       ierr = VecPointwiseDivide(shell->right_work,D,shell->right);CHKERRQ(ierr);
428     }
429     if (!shell->dshift) {
430       ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->dshift);CHKERRQ(ierr);
431       ierr = VecCopy(shell->dshift,shell->right_work);CHKERRQ(ierr);
432     } else {
433       ierr = VecAXPY(shell->dshift,1.0,shell->right_work);CHKERRQ(ierr);
434     }
435   } else {
436     ierr = VecAXPY(shell->dshift,1.0,D);CHKERRQ(ierr);
437   }
438   PetscFunctionReturn(0);
439 }
440 
441 PetscErrorCode MatScale_Shell(Mat Y,PetscScalar a)
442 {
443   Mat_Shell      *shell = (Mat_Shell*)Y->data;
444   PetscErrorCode ierr;
445 
446   PetscFunctionBegin;
447   shell->vscale *= a;
448   shell->vshift *= a;
449   if (shell->dshift) {
450     ierr = VecScale(shell->dshift,a);CHKERRQ(ierr);
451   }
452   shell->axpy_vscale *= a;
453   PetscFunctionReturn(0);
454 }
455 
456 static PetscErrorCode MatDiagonalScale_Shell(Mat Y,Vec left,Vec right)
457 {
458   Mat_Shell      *shell = (Mat_Shell*)Y->data;
459   PetscErrorCode ierr;
460 
461   PetscFunctionBegin;
462   if (shell->axpy) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_SUP,"Cannot diagonal scale MATSHELL after MatAXPY operation");
463   if (left) {
464     if (!shell->left) {
465       ierr = VecDuplicate(left,&shell->left);CHKERRQ(ierr);
466       ierr = VecCopy(left,shell->left);CHKERRQ(ierr);
467     } else {
468       ierr = VecPointwiseMult(shell->left,shell->left,left);CHKERRQ(ierr);
469     }
470   }
471   if (right) {
472     if (!shell->right) {
473       ierr = VecDuplicate(right,&shell->right);CHKERRQ(ierr);
474       ierr = VecCopy(right,shell->right);CHKERRQ(ierr);
475     } else {
476       ierr = VecPointwiseMult(shell->right,shell->right,right);CHKERRQ(ierr);
477     }
478   }
479   PetscFunctionReturn(0);
480 }
481 
482 PetscErrorCode MatAssemblyEnd_Shell(Mat Y,MatAssemblyType t)
483 {
484   Mat_Shell      *shell = (Mat_Shell*)Y->data;
485   PetscErrorCode ierr;
486 
487   PetscFunctionBegin;
488   if (t == MAT_FINAL_ASSEMBLY) {
489     shell->vshift = 0.0;
490     shell->vscale = 1.0;
491     ierr = VecDestroy(&shell->dshift);CHKERRQ(ierr);
492     ierr = VecDestroy(&shell->left);CHKERRQ(ierr);
493     ierr = VecDestroy(&shell->right);CHKERRQ(ierr);
494     ierr = MatDestroy(&shell->axpy);CHKERRQ(ierr);
495   }
496   PetscFunctionReturn(0);
497 }
498 
499 PETSC_INTERN PetscErrorCode MatConvert_Shell(Mat, MatType,MatReuse,Mat*);
500 
501 static PetscErrorCode MatMissingDiagonal_Shell(Mat A,PetscBool  *missing,PetscInt *d)
502 {
503   PetscFunctionBegin;
504   *missing = PETSC_FALSE;
505   PetscFunctionReturn(0);
506 }
507 
508 PetscErrorCode MatAXPY_Shell(Mat Y,PetscScalar a,Mat X,MatStructure str)
509 {
510   Mat_Shell      *shell = (Mat_Shell*)Y->data;
511   PetscErrorCode ierr;
512 
513   PetscFunctionBegin;
514   ierr = PetscObjectReference((PetscObject)X);CHKERRQ(ierr);
515   ierr = MatDestroy(&shell->axpy);CHKERRQ(ierr);
516   shell->axpy        = X;
517   shell->axpy_vscale = a;
518   PetscFunctionReturn(0);
519 }
520 
521 static struct _MatOps MatOps_Values = {0,
522                                        0,
523                                        0,
524                                        0,
525                                 /* 4*/ MatMultAdd_Shell,
526                                        0,
527                                        MatMultTransposeAdd_Shell,
528                                        0,
529                                        0,
530                                        0,
531                                 /*10*/ 0,
532                                        0,
533                                        0,
534                                        0,
535                                        0,
536                                 /*15*/ 0,
537                                        0,
538                                        MatGetDiagonal_Shell,
539                                        MatDiagonalScale_Shell,
540                                        0,
541                                 /*20*/ 0,
542                                        MatAssemblyEnd_Shell,
543                                        0,
544                                        0,
545                                 /*24*/ 0,
546                                        0,
547                                        0,
548                                        0,
549                                        0,
550                                 /*29*/ 0,
551                                        0,
552                                        0,
553                                        0,
554                                        0,
555                                 /*34*/ MatDuplicate_Shell,
556                                        0,
557                                        0,
558                                        0,
559                                        0,
560                                 /*39*/ MatAXPY_Shell,
561                                        0,
562                                        0,
563                                        0,
564                                        MatCopy_Shell,
565                                 /*44*/ 0,
566                                        MatScale_Shell,
567                                        MatShift_Shell,
568                                        MatDiagonalSet_Shell,
569                                        0,
570                                 /*49*/ 0,
571                                        0,
572                                        0,
573                                        0,
574                                        0,
575                                 /*54*/ 0,
576                                        0,
577                                        0,
578                                        0,
579                                        0,
580                                 /*59*/ 0,
581                                        MatDestroy_Shell,
582                                        0,
583                                        0,
584                                        0,
585                                 /*64*/ 0,
586                                        0,
587                                        0,
588                                        0,
589                                        0,
590                                 /*69*/ 0,
591                                        0,
592                                        MatConvert_Shell,
593                                        0,
594                                        0,
595                                 /*74*/ 0,
596                                        0,
597                                        0,
598                                        0,
599                                        0,
600                                 /*79*/ 0,
601                                        0,
602                                        0,
603                                        0,
604                                        0,
605                                 /*84*/ 0,
606                                        0,
607                                        0,
608                                        0,
609                                        0,
610                                 /*89*/ 0,
611                                        0,
612                                        0,
613                                        0,
614                                        0,
615                                 /*94*/ 0,
616                                        0,
617                                        0,
618                                        0,
619                                        0,
620                                 /*99*/ 0,
621                                        0,
622                                        0,
623                                        0,
624                                        0,
625                                /*104*/ 0,
626                                        0,
627                                        0,
628                                        0,
629                                        0,
630                                /*109*/ 0,
631                                        0,
632                                        0,
633                                        0,
634                                        MatMissingDiagonal_Shell,
635                                /*114*/ 0,
636                                        0,
637                                        0,
638                                        0,
639                                        0,
640                                /*119*/ 0,
641                                        0,
642                                        0,
643                                        0,
644                                        0,
645                                /*124*/ 0,
646                                        0,
647                                        0,
648                                        0,
649                                        0,
650                                /*129*/ 0,
651                                        0,
652                                        0,
653                                        0,
654                                        0,
655                                /*134*/ 0,
656                                        0,
657                                        0,
658                                        0,
659                                        0,
660                                /*139*/ 0,
661                                        0,
662                                        0
663 };
664 
665 /*MC
666    MATSHELL - MATSHELL = "shell" - A matrix type to be used to define your own matrix type -- perhaps matrix free.
667 
668   Level: advanced
669 
670 .seealso: MatCreateShell()
671 M*/
672 
673 PETSC_EXTERN PetscErrorCode MatCreate_Shell(Mat A)
674 {
675   Mat_Shell      *b;
676   PetscErrorCode ierr;
677 
678   PetscFunctionBegin;
679   ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
680 
681   ierr    = PetscNewLog(A,&b);CHKERRQ(ierr);
682   A->data = (void*)b;
683 
684   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
685   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
686 
687   b->ctx                 = 0;
688   b->vshift              = 0.0;
689   b->vscale              = 1.0;
690   b->managescalingshifts = PETSC_TRUE;
691   A->assembled           = PETSC_TRUE;
692   A->preallocated        = PETSC_FALSE;
693 
694   ierr = PetscObjectChangeTypeName((PetscObject)A,MATSHELL);CHKERRQ(ierr);
695   PetscFunctionReturn(0);
696 }
697 
698 /*@C
699    MatCreateShell - Creates a new matrix class for use with a user-defined
700    private data storage format.
701 
702   Collective on MPI_Comm
703 
704    Input Parameters:
705 +  comm - MPI communicator
706 .  m - number of local rows (must be given)
707 .  n - number of local columns (must be given)
708 .  M - number of global rows (may be PETSC_DETERMINE)
709 .  N - number of global columns (may be PETSC_DETERMINE)
710 -  ctx - pointer to data needed by the shell matrix routines
711 
712    Output Parameter:
713 .  A - the matrix
714 
715    Level: advanced
716 
717   Usage:
718 $    extern int mult(Mat,Vec,Vec);
719 $    MatCreateShell(comm,m,n,M,N,ctx,&mat);
720 $    MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult);
721 $    [ Use matrix for operations that have been set ]
722 $    MatDestroy(mat);
723 
724    Notes:
725    The shell matrix type is intended to provide a simple class to use
726    with KSP (such as, for use with matrix-free methods). You should not
727    use the shell type if you plan to define a complete matrix class.
728 
729    Fortran Notes: To use this from Fortran with a ctx you must write an interface definition for this
730     function and for MatShellGetContext() that tells Fortran the Fortran derived data type you are passing
731     in as the ctx argument.
732 
733    PETSc requires that matrices and vectors being used for certain
734    operations are partitioned accordingly.  For example, when
735    creating a shell matrix, A, that supports parallel matrix-vector
736    products using MatMult(A,x,y) the user should set the number
737    of local matrix rows to be the number of local elements of the
738    corresponding result vector, y. Note that this is information is
739    required for use of the matrix interface routines, even though
740    the shell matrix may not actually be physically partitioned.
741    For example,
742 
743    MATSHELL handles MatShift(), MatDiagonalSet(), MatDiagonalScale(), and MatScale() internally so these
744    operations cannot be overwritten unless MatShellSetManageScalingShifts() is called.
745 
746 $
747 $     Vec x, y
748 $     extern int mult(Mat,Vec,Vec);
749 $     Mat A
750 $
751 $     VecCreateMPI(comm,PETSC_DECIDE,M,&y);
752 $     VecCreateMPI(comm,PETSC_DECIDE,N,&x);
753 $     VecGetLocalSize(y,&m);
754 $     VecGetLocalSize(x,&n);
755 $     MatCreateShell(comm,m,n,M,N,ctx,&A);
756 $     MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult);
757 $     MatMult(A,x,y);
758 $     MatDestroy(A);
759 $     VecDestroy(y); VecDestroy(x);
760 $
761 
762     For rectangular matrices do all the scalings and shifts make sense?
763 
764     Developers Notes: Regarding shifting and scaling. The general form is
765 
766           diag(left)(vscale*A + diag(dshift) + vshift I)diag(right)
767 
768       The order you apply the operations is important. For example if you have a dshift then
769       apply a MatScale(s) you get s*vscale*A + s*diag(shift). But if you first scale and then shift
770       you get s*vscale*A + diag(shift)
771 
772           A is the user provided function.
773 
774 .keywords: matrix, shell, create
775 
776 .seealso: MatShellSetOperation(), MatHasOperation(), MatShellGetContext(), MatShellSetContext(), MATSHELL, MatShellSetManageScalingShifts()
777 @*/
778 PetscErrorCode  MatCreateShell(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx,Mat *A)
779 {
780   PetscErrorCode ierr;
781 
782   PetscFunctionBegin;
783   ierr = MatCreate(comm,A);CHKERRQ(ierr);
784   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
785   ierr = MatSetType(*A,MATSHELL);CHKERRQ(ierr);
786   ierr = MatShellSetContext(*A,ctx);CHKERRQ(ierr);
787   ierr = MatSetUp(*A);CHKERRQ(ierr);
788   PetscFunctionReturn(0);
789 }
790 
791 /*@
792     MatShellSetContext - sets the context for a shell matrix
793 
794    Logically Collective on Mat
795 
796     Input Parameters:
797 +   mat - the shell matrix
798 -   ctx - the context
799 
800    Level: advanced
801 
802    Fortran Notes: To use this from Fortran you must write a Fortran interface definition for this
803     function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.
804 
805 .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation()
806 @*/
807 PetscErrorCode  MatShellSetContext(Mat mat,void *ctx)
808 {
809   Mat_Shell      *shell = (Mat_Shell*)mat->data;
810   PetscErrorCode ierr;
811   PetscBool      flg;
812 
813   PetscFunctionBegin;
814   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
815   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
816   if (flg) {
817     shell->ctx = ctx;
818   } else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot attach context to non-shell matrix");
819   PetscFunctionReturn(0);
820 }
821 
822 /*@
823     MatShellSetManageScalingShifts - Allows the user to control the scaling and shift operations of the MATSHELL. Must be called immediately
824           after MatCreateShell()
825 
826    Logically Collective on Mat
827 
828     Input Parameter:
829 .   mat - the shell matrix
830 
831   Level: advanced
832 
833 .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext(), MatShellSetOperation()
834 @*/
835 PetscErrorCode MatShellSetManageScalingShifts(Mat A)
836 {
837   PetscErrorCode ierr;
838   Mat_Shell      *shell = (Mat_Shell*)A->data;
839   PetscBool      flg;
840 
841   PetscFunctionBegin;
842   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
843   ierr = PetscObjectTypeCompare((PetscObject)A,MATSHELL,&flg);CHKERRQ(ierr);
844   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Can only use with MATSHELL matrices");
845   shell->managescalingshifts = PETSC_FALSE;
846   A->ops->diagonalscale = 0;
847   A->ops->scale         = 0;
848   A->ops->shift         = 0;
849   A->ops->diagonalset   = 0;
850   PetscFunctionReturn(0);
851 }
852 
853 /*@C
854     MatShellTestMult - Compares the multiply routine provided to the MATSHELL with differencing on a given function.
855 
856    Logically Collective on Mat
857 
858     Input Parameters:
859 +   mat - the shell matrix
860 .   f - the function
861 .   base - differences are computed around this vector, see MatMFFDSetBase(), for Jacobians this is the point at which the Jacobian is being evaluated
862 -   ctx - an optional context for the function
863 
864    Output Parameter:
865 .   flg - PETSC_TRUE if the multiply is likely correct
866 
867    Options Database:
868 .   -mat_shell_test_mult_view - print if any differences are detected between the products and print the difference
869 
870    Level: advanced
871 
872    Fortran Notes: Not supported from Fortran
873 
874 .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellTestMultTranspose()
875 @*/
876 PetscErrorCode  MatShellTestMult(Mat mat,PetscErrorCode (*f)(void*,Vec,Vec),Vec base,void *ctx,PetscBool *flg)
877 {
878   PetscErrorCode ierr;
879   PetscInt       m,n;
880   Mat            mf,Dmf,Dmat,Ddiff;
881   PetscReal      Diffnorm,Dmfnorm;
882   PetscBool      v = PETSC_FALSE, flag = PETSC_TRUE;
883 
884   PetscFunctionBegin;
885   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
886   ierr = PetscOptionsHasName(NULL,((PetscObject)mat)->prefix,"-mat_shell_test_mult_view",&v);CHKERRQ(ierr);
887   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
888   ierr = MatCreateMFFD(PetscObjectComm((PetscObject)mat),m,n,PETSC_DECIDE,PETSC_DECIDE,&mf);CHKERRQ(ierr);
889   ierr = MatMFFDSetFunction(mf,f,ctx);CHKERRQ(ierr);
890   ierr = MatMFFDSetBase(mf,base,NULL);CHKERRQ(ierr);
891 
892   ierr = MatComputeExplicitOperator(mf,&Dmf);CHKERRQ(ierr);
893   ierr = MatComputeExplicitOperator(mat,&Dmat);CHKERRQ(ierr);
894 
895   ierr = MatDuplicate(Dmat,MAT_COPY_VALUES,&Ddiff);CHKERRQ(ierr);
896   ierr = MatAXPY(Ddiff,-1.0,Dmf,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
897   ierr = MatNorm(Ddiff,NORM_FROBENIUS,&Diffnorm);CHKERRQ(ierr);
898   ierr = MatNorm(Dmf,NORM_FROBENIUS,&Dmfnorm);CHKERRQ(ierr);
899   if (Diffnorm/Dmfnorm > 10*PETSC_SQRT_MACHINE_EPSILON) {
900     flag = PETSC_FALSE;
901     if (v) {
902       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));
903       ierr = MatViewFromOptions(Ddiff,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr);
904       ierr = MatViewFromOptions(Dmf,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr);
905       ierr = MatViewFromOptions(Dmat,(PetscObject)mat,"-mat_shell_test_mult_view");CHKERRQ(ierr);
906     }
907   } else if (v) {
908     ierr = PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL and matrix free multiple appear to produce the same results\n");
909   }
910   if (flg) *flg = flag;
911   ierr = MatDestroy(&Ddiff);CHKERRQ(ierr);
912   ierr = MatDestroy(&mf);CHKERRQ(ierr);
913   ierr = MatDestroy(&Dmf);CHKERRQ(ierr);
914   ierr = MatDestroy(&Dmat);CHKERRQ(ierr);
915   PetscFunctionReturn(0);
916 }
917 
918 /*@C
919     MatShellTestMultTranpose - Compares the multiply transpose routine provided to the MATSHELL with differencing on a given function.
920 
921    Logically Collective on Mat
922 
923     Input Parameters:
924 +   mat - the shell matrix
925 .   f - the function
926 .   base - differences are computed around this vector, see MatMFFDSetBase(), for Jacobians this is the point at which the Jacobian is being evaluated
927 -   ctx - an optional context for the function
928 
929    Output Parameter:
930 .   flg - PETSC_TRUE if the multiply is likely correct
931 
932    Options Database:
933 .   -mat_shell_test_mult_view - print if any differences are detected between the products and print the difference
934 
935    Level: advanced
936 
937    Fortran Notes: Not supported from Fortran
938 
939 .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellTestMult()
940 @*/
941 PetscErrorCode  MatShellTestMultTranspose(Mat mat,PetscErrorCode (*f)(void*,Vec,Vec),Vec base,void *ctx,PetscBool *flg)
942 {
943   PetscErrorCode ierr;
944   Vec            x,y,z;
945   PetscInt       m,n,M,N;
946   Mat            mf,Dmf,Dmat,Ddiff;
947   PetscReal      Diffnorm,Dmfnorm;
948   PetscBool      v = PETSC_FALSE, flag = PETSC_TRUE;
949 
950   PetscFunctionBegin;
951   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
952   ierr = PetscOptionsHasName(NULL,((PetscObject)mat)->prefix,"-mat_shell_test_mult_transpose_view",&v);CHKERRQ(ierr);
953   ierr = MatCreateVecs(mat,&x,&y);CHKERRQ(ierr);
954   ierr = VecDuplicate(y,&z);CHKERRQ(ierr);
955   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
956   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
957   ierr = MatCreateMFFD(PetscObjectComm((PetscObject)mat),m,n,M,N,&mf);CHKERRQ(ierr);
958   ierr = MatMFFDSetFunction(mf,f,ctx);CHKERRQ(ierr);
959   ierr = MatMFFDSetBase(mf,base,NULL);CHKERRQ(ierr);
960   ierr = MatComputeExplicitOperator(mf,&Dmf);CHKERRQ(ierr);
961   ierr = MatTranspose(Dmf,MAT_INPLACE_MATRIX,&Dmf);CHKERRQ(ierr);
962   ierr = MatComputeExplicitOperatorTranspose(mat,&Dmat);CHKERRQ(ierr);
963 
964   ierr = MatDuplicate(Dmat,MAT_COPY_VALUES,&Ddiff);CHKERRQ(ierr);
965   ierr = MatAXPY(Ddiff,-1.0,Dmf,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
966   ierr = MatNorm(Ddiff,NORM_FROBENIUS,&Diffnorm);CHKERRQ(ierr);
967   ierr = MatNorm(Dmf,NORM_FROBENIUS,&Dmfnorm);CHKERRQ(ierr);
968   if (Diffnorm/Dmfnorm > 10*PETSC_SQRT_MACHINE_EPSILON) {
969     flag = PETSC_FALSE;
970     if (v) {
971       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));
972       ierr = MatViewFromOptions(Ddiff,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr);
973       ierr = MatViewFromOptions(Dmf,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr);
974       ierr = MatViewFromOptions(Dmat,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");CHKERRQ(ierr);
975     }
976   } else if (v) {
977     ierr = PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL transpose and matrix free multiple appear to produce the same results\n");
978   }
979   if (flg) *flg = flag;
980   ierr = MatDestroy(&mf);CHKERRQ(ierr);
981   ierr = MatDestroy(&Dmat);CHKERRQ(ierr);
982   ierr = MatDestroy(&Ddiff);CHKERRQ(ierr);
983   ierr = MatDestroy(&Dmf);CHKERRQ(ierr);
984   ierr = VecDestroy(&x);CHKERRQ(ierr);
985   ierr = VecDestroy(&y);CHKERRQ(ierr);
986   ierr = VecDestroy(&z);CHKERRQ(ierr);
987   PetscFunctionReturn(0);
988 }
989 
990 /*@C
991     MatShellSetOperation - Allows user to set a matrix operation for
992                            a shell matrix.
993     MatShellSetOperation - Allows user to set a matrix operation for a shell matrix.
994 
995    Logically Collective on Mat
996 
997     Input Parameters:
998 +   mat - the shell matrix
999 .   op - the name of the operation
1000 -   f - the function that provides the operation.
1001 
1002    Level: advanced
1003 
1004     Usage:
1005 $      extern PetscErrorCode usermult(Mat,Vec,Vec);
1006 $      ierr = MatCreateShell(comm,m,n,M,N,ctx,&A);
1007 $      ierr = MatShellSetOperation(A,MATOP_MULT,(void(*)(void))usermult);
1008 
1009     Notes:
1010     See the file include/petscmat.h for a complete list of matrix
1011     operations, which all have the form MATOP_<OPERATION>, where
1012     <OPERATION> is the name (in all capital letters) of the
1013     user interface routine (e.g., MatMult() -> MATOP_MULT).
1014 
1015     All user-provided functions (except for MATOP_DESTROY) should have the same calling
1016     sequence as the usual matrix interface routines, since they
1017     are intended to be accessed via the usual matrix interface
1018     routines, e.g.,
1019 $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)
1020 
1021     In particular each function MUST return an error code of 0 on success and
1022     nonzero on failure.
1023 
1024     Within each user-defined routine, the user should call
1025     MatShellGetContext() to obtain the user-defined context that was
1026     set by MatCreateShell().
1027 
1028     Fortran Notes: For MatCreateVecs() the user code should check if the input left or right matrix is -1 and in that case not
1029        generate a matrix. See src/mat/examples/tests/ex120f.F
1030 
1031     Use MatSetOperation() to set an operation for any matrix type
1032 
1033 .keywords: matrix, shell, set, operation
1034 
1035 .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext(), MatSetOperation(), MatShellSetManageScalingShifts()
1036 @*/
1037 PetscErrorCode  MatShellSetOperation(Mat mat,MatOperation op,void (*f)(void))
1038 {
1039   PetscErrorCode ierr;
1040   PetscBool      flg;
1041   Mat_Shell      *shell = (Mat_Shell*)mat->data;
1042 
1043   PetscFunctionBegin;
1044   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
1045   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
1046   if (!flg) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Can only use with MATSHELL matrices");
1047   switch (op) {
1048   case MATOP_DESTROY:
1049     shell->ops->destroy = (PetscErrorCode (*)(Mat))f;
1050     break;
1051   case MATOP_DIAGONAL_SET:
1052   case MATOP_DIAGONAL_SCALE:
1053   case MATOP_SHIFT:
1054   case MATOP_SCALE:
1055   case MATOP_AXPY:
1056     if (shell->managescalingshifts) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"MATSHELL is managing scalings and shifts, see MatShellSetManageScalingShifts()");
1057     (((void(**)(void))mat->ops)[op]) = f;
1058     break;
1059   case MATOP_GET_DIAGONAL:
1060     if (shell->managescalingshifts) shell->ops->getdiagonal = (PetscErrorCode (*)(Mat,Vec))f;
1061     else mat->ops->getdiagonal = (PetscErrorCode (*)(Mat,Vec))f;
1062   case MATOP_VIEW:
1063     if (!mat->ops->viewnative) {
1064       mat->ops->viewnative = mat->ops->view;
1065     }
1066     mat->ops->view = (PetscErrorCode (*)(Mat,PetscViewer))f;
1067     break;
1068   case MATOP_MULT:
1069     if (shell->managescalingshifts){
1070       shell->ops->mult = (PetscErrorCode (*)(Mat,Vec,Vec))f;
1071       mat->ops->mult   = MatMult_Shell;
1072     } else mat->ops->mult = (PetscErrorCode (*)(Mat,Vec,Vec))f;
1073     break;
1074   case MATOP_MULT_TRANSPOSE:
1075     if (shell->managescalingshifts) {
1076       shell->ops->multtranspose = (PetscErrorCode (*)(Mat,Vec,Vec))f;
1077       mat->ops->multtranspose   = MatMultTranspose_Shell;
1078     } else mat->ops->multtranspose = (PetscErrorCode (*)(Mat,Vec,Vec))f;
1079     break;
1080   default:
1081     (((void(**)(void))mat->ops)[op]) = f;
1082   }
1083   PetscFunctionReturn(0);
1084 }
1085 
1086 /*@C
1087     MatShellGetOperation - Gets a matrix function for a shell matrix.
1088 
1089     Not Collective
1090 
1091     Input Parameters:
1092 +   mat - the shell matrix
1093 -   op - the name of the operation
1094 
1095     Output Parameter:
1096 .   f - the function that provides the operation.
1097 
1098     Level: advanced
1099 
1100     Notes:
1101     See the file include/petscmat.h for a complete list of matrix
1102     operations, which all have the form MATOP_<OPERATION>, where
1103     <OPERATION> is the name (in all capital letters) of the
1104     user interface routine (e.g., MatMult() -> MATOP_MULT).
1105 
1106     All user-provided functions have the same calling
1107     sequence as the usual matrix interface routines, since they
1108     are intended to be accessed via the usual matrix interface
1109     routines, e.g.,
1110 $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)
1111 
1112     Within each user-defined routine, the user should call
1113     MatShellGetContext() to obtain the user-defined context that was
1114     set by MatCreateShell().
1115 
1116 .keywords: matrix, shell, set, operation
1117 
1118 .seealso: MatCreateShell(), MatShellGetContext(), MatShellSetOperation(), MatShellSetContext()
1119 @*/
1120 PetscErrorCode  MatShellGetOperation(Mat mat,MatOperation op,void(**f)(void))
1121 {
1122   PetscErrorCode ierr;
1123   PetscBool      flg;
1124 
1125   PetscFunctionBegin;
1126   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
1127   switch (op) {
1128   case MATOP_DESTROY:
1129     ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
1130     if (flg) {
1131       Mat_Shell *shell = (Mat_Shell*)mat->data;
1132       *f = (void (*)(void))shell->ops->destroy;
1133     } else {
1134       *f = (void (*)(void))mat->ops->destroy;
1135     }
1136     break;
1137   case MATOP_DIAGONAL_SET:
1138     ierr = PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
1139     if (flg) {
1140       Mat_Shell *shell = (Mat_Shell*)mat->data;
1141       *f = (void (*)(void))shell->ops->diagonalset;
1142     } else {
1143       *f = (void (*)(void))mat->ops->diagonalset;
1144     }
1145     break;
1146   case MATOP_VIEW:
1147     *f = (void (*)(void))mat->ops->view;
1148     break;
1149   default:
1150     *f = (((void (**)(void))mat->ops)[op]);
1151   }
1152   PetscFunctionReturn(0);
1153 }
1154 
1155 /*@C
1156     MatSetOperation - Allows user to set a matrix operation for any matrix type
1157 
1158    Logically Collective on Mat
1159 
1160     Input Parameters:
1161 +   mat - the shell matrix
1162 .   op - the name of the operation
1163 -   f - the function that provides the operation.
1164 
1165    Level: developer
1166 
1167     Usage:
1168 $      extern PetscErrorCode usermult(Mat,Vec,Vec);
1169 $      ierr = MatCreateXXX(comm,...&A);
1170 $      ierr = MatShellSetOperation(A,MATOP_MULT,(void(*)(void))usermult);
1171 
1172     Notes:
1173     See the file include/petscmat.h for a complete list of matrix
1174     operations, which all have the form MATOP_<OPERATION>, where
1175     <OPERATION> is the name (in all capital letters) of the
1176     user interface routine (e.g., MatMult() -> MATOP_MULT).
1177 
1178     All user-provided functions (except for MATOP_DESTROY) should have the same calling
1179     sequence as the usual matrix interface routines, since they
1180     are intended to be accessed via the usual matrix interface
1181     routines, e.g.,
1182 $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)
1183 
1184     In particular each function MUST return an error code of 0 on success and
1185     nonzero on failure.
1186 
1187     This routine is distinct from MatShellSetOperation() in that it can be called on any matrix type.
1188 
1189 .keywords: matrix, shell, set, operation
1190 
1191 .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext(), MatShellSetOperation()
1192 @*/
1193 PetscErrorCode  MatSetOperation(Mat mat,MatOperation op,void (*f)(void))
1194 {
1195   PetscFunctionBegin;
1196   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
1197   (((void(**)(void))mat->ops)[op]) = f;
1198   PetscFunctionReturn(0);
1199 }
1200