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