xref: /petsc/src/mat/impls/shell/shell.c (revision f73d5cc4fab9985badaa1a29b152f456bf59ad34)
1 
2 /*
3    This provides a simple shell for Fortran (and C programmers) to
4   create a very simple matrix class for use with KSP without coding
5   much of anything.
6 */
7 
8 #include <petsc-private/matimpl.h>        /*I "petscmat.h" I*/
9 #include <petsc-private/vecimpl.h>
10 
11 typedef struct {
12   PetscErrorCode (*destroy)(Mat);
13   PetscErrorCode (*mult)(Mat,Vec,Vec);
14   PetscErrorCode (*multtranspose)(Mat,Vec,Vec);
15   PetscErrorCode (*getdiagonal)(Mat,Vec);
16   PetscScalar    vscale,vshift;
17   Vec            dshift;
18   Vec            left,right;
19   Vec            dshift_owned,left_owned,right_owned;
20   Vec            left_work,right_work;
21   PetscBool      usingscaled;
22   void           *ctx;
23 } Mat_Shell;
24 /*
25  The most general expression for the matrix is
26 
27  S = L (a A + B) R
28 
29  where
30  A is the matrix defined by the user's function
31  a is a scalar multiple
32  L is left scaling
33  R is right scaling
34  B is a diagonal shift defined by
35    diag(dshift) if the vector dshift is non-NULL
36    vscale*identity otherwise
37 
38  The following identities apply:
39 
40  Scale by c:
41   c [L (a A + B) R] = L [(a c) A + c B] R
42 
43  Shift by c:
44   [L (a A + B) R] + c = L [a A + (B + c Linv Rinv)] R
45 
46  Diagonal scaling is achieved by simply multiplying with existing L and R vectors
47 
48  In the data structure:
49 
50  vscale=1.0  means no special scaling will be applied
51  dshift=NULL means a constant diagonal shift (fall back to vshift)
52  vshift=0.0  means no constant diagonal shift, note that vshift is only used if dshift is NULL
53 */
54 
55 static PetscErrorCode MatMult_Shell(Mat,Vec,Vec);
56 static PetscErrorCode MatMultTranspose_Shell(Mat,Vec,Vec);
57 static PetscErrorCode MatGetDiagonal_Shell(Mat,Vec);
58 
59 #undef __FUNCT__
60 #define __FUNCT__ "MatShellUseScaledMethods"
61 static PetscErrorCode MatShellUseScaledMethods(Mat Y)
62 {
63   Mat_Shell *shell = (Mat_Shell*)Y->data;
64 
65   PetscFunctionBegin;
66   if (shell->usingscaled) PetscFunctionReturn(0);
67   shell->mult  = Y->ops->mult;
68   Y->ops->mult = MatMult_Shell;
69   if (Y->ops->multtranspose) {
70     shell->multtranspose  = Y->ops->multtranspose;
71     Y->ops->multtranspose = MatMultTranspose_Shell;
72   }
73   if (Y->ops->getdiagonal) {
74     shell->getdiagonal  = Y->ops->getdiagonal;
75     Y->ops->getdiagonal = MatGetDiagonal_Shell;
76   }
77   shell->usingscaled = PETSC_TRUE;
78   PetscFunctionReturn(0);
79 }
80 
81 #undef __FUNCT__
82 #define __FUNCT__ "MatShellPreScaleLeft"
83 static PetscErrorCode MatShellPreScaleLeft(Mat A,Vec x,Vec *xx)
84 {
85   Mat_Shell *shell = (Mat_Shell*)A->data;
86   PetscErrorCode ierr;
87 
88   PetscFunctionBegin;
89   *xx = PETSC_NULL;
90   if (!shell->left) {
91     *xx = x;
92   } else {
93     if (!shell->left_work) {ierr = VecDuplicate(shell->left,&shell->left_work);CHKERRQ(ierr);}
94     ierr = VecPointwiseMult(shell->left_work,x,shell->left);CHKERRQ(ierr);
95     *xx = shell->left_work;
96   }
97   PetscFunctionReturn(0);
98 }
99 
100 #undef __FUNCT__
101 #define __FUNCT__ "MatShellPreScaleRight"
102 static PetscErrorCode MatShellPreScaleRight(Mat A,Vec x,Vec *xx)
103 {
104   Mat_Shell *shell = (Mat_Shell*)A->data;
105   PetscErrorCode ierr;
106 
107   PetscFunctionBegin;
108   *xx = PETSC_NULL;
109   if (!shell->right) {
110     *xx = x;
111   } else {
112     if (!shell->right_work) {ierr = VecDuplicate(shell->right,&shell->right_work);CHKERRQ(ierr);}
113     ierr = VecPointwiseMult(shell->right_work,x,shell->right);CHKERRQ(ierr);
114     *xx = shell->right_work;
115   }
116   PetscFunctionReturn(0);
117 }
118 
119 #undef __FUNCT__
120 #define __FUNCT__ "MatShellPostScaleLeft"
121 static PetscErrorCode MatShellPostScaleLeft(Mat A,Vec x)
122 {
123   Mat_Shell *shell = (Mat_Shell*)A->data;
124   PetscErrorCode ierr;
125 
126   PetscFunctionBegin;
127   if (shell->left) {ierr = VecPointwiseMult(x,x,shell->left);CHKERRQ(ierr);}
128   PetscFunctionReturn(0);
129 }
130 
131 #undef __FUNCT__
132 #define __FUNCT__ "MatShellPostScaleRight"
133 static PetscErrorCode MatShellPostScaleRight(Mat A,Vec x)
134 {
135   Mat_Shell *shell = (Mat_Shell*)A->data;
136   PetscErrorCode ierr;
137 
138   PetscFunctionBegin;
139   if (shell->right) {ierr = VecPointwiseMult(x,x,shell->right);CHKERRQ(ierr);}
140   PetscFunctionReturn(0);
141 }
142 
143 #undef __FUNCT__
144 #define __FUNCT__ "MatShellShiftAndScale"
145 static PetscErrorCode MatShellShiftAndScale(Mat A,Vec X,Vec Y)
146 {
147   Mat_Shell *shell = (Mat_Shell*)A->data;
148   PetscErrorCode ierr;
149 
150   PetscFunctionBegin;
151   if (shell->dshift) {          /* get arrays because there is no VecPointwiseMultAdd() */
152     PetscInt          i,m;
153     const PetscScalar *x,*d;
154     PetscScalar       *y;
155     ierr = VecGetLocalSize(X,&m);CHKERRQ(ierr);
156     ierr = VecGetArrayRead(shell->dshift,&d);CHKERRQ(ierr);
157     ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
158     ierr = VecGetArray(Y,&y);CHKERRQ(ierr);
159     for (i=0; i<m; i++) y[i] = shell->vscale*y[i] + d[i]*x[i];
160     ierr = VecRestoreArrayRead(shell->dshift,&d);CHKERRQ(ierr);
161     ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
162     ierr = VecRestoreArray(Y,&y);CHKERRQ(ierr);
163   } else if (PetscAbsScalar(shell->vshift) != 0) {
164     ierr = VecAXPBY(Y,shell->vshift,shell->vscale,X);CHKERRQ(ierr);
165   } else if (shell->vscale != 1.0) {
166     ierr = VecScale(Y,shell->vscale);CHKERRQ(ierr);
167   }
168   PetscFunctionReturn(0);
169 }
170 
171 #undef __FUNCT__
172 #define __FUNCT__ "MatShellGetContext"
173 /*@
174     MatShellGetContext - Returns the user-provided context associated with a shell matrix.
175 
176     Not Collective
177 
178     Input Parameter:
179 .   mat - the matrix, should have been created with MatCreateShell()
180 
181     Output Parameter:
182 .   ctx - the user provided context
183 
184     Level: advanced
185 
186     Notes:
187     This routine is intended for use within various shell matrix routines,
188     as set with MatShellSetOperation().
189 
190 .keywords: matrix, shell, get, context
191 
192 .seealso: MatCreateShell(), MatShellSetOperation(), MatShellSetContext()
193 @*/
194 PetscErrorCode  MatShellGetContext(Mat mat,void *ctx)
195 {
196   PetscErrorCode ierr;
197   PetscBool      flg;
198 
199   PetscFunctionBegin;
200   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
201   PetscValidPointer(ctx,2);
202   ierr = PetscTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
203   if (flg) *(void**)ctx = ((Mat_Shell*)(mat->data))->ctx;
204   else SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_SUP,"Cannot get context from non-shell matrix");
205   PetscFunctionReturn(0);
206 }
207 
208 #undef __FUNCT__
209 #define __FUNCT__ "MatDestroy_Shell"
210 PetscErrorCode MatDestroy_Shell(Mat mat)
211 {
212   PetscErrorCode ierr;
213   Mat_Shell      *shell = (Mat_Shell*)mat->data;
214 
215   PetscFunctionBegin;
216   if (shell->destroy) {
217     ierr = (*shell->destroy)(mat);CHKERRQ(ierr);
218   }
219   ierr = VecDestroy(&shell->left_owned);CHKERRQ(ierr);
220   ierr = VecDestroy(&shell->right_owned);CHKERRQ(ierr);
221   ierr = VecDestroy(&shell->dshift_owned);CHKERRQ(ierr);
222   ierr = VecDestroy(&shell->left_work);CHKERRQ(ierr);
223   ierr = VecDestroy(&shell->right_work);CHKERRQ(ierr);
224   ierr = PetscFree(mat->data);CHKERRQ(ierr);
225   PetscFunctionReturn(0);
226 }
227 
228 #undef __FUNCT__
229 #define __FUNCT__ "MatMult_Shell"
230 PetscErrorCode MatMult_Shell(Mat A,Vec x,Vec y)
231 {
232   Mat_Shell      *shell = (Mat_Shell*)A->data;
233   PetscErrorCode ierr;
234   Vec            xx;
235 
236   PetscFunctionBegin;
237   ierr = MatShellPreScaleRight(A,x,&xx);CHKERRQ(ierr);
238   ierr = (*shell->mult)(A,xx,y);CHKERRQ(ierr);
239   ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr);
240   ierr = MatShellPostScaleLeft(A,y);CHKERRQ(ierr);
241   PetscFunctionReturn(0);
242 }
243 
244 #undef __FUNCT__
245 #define __FUNCT__ "MatMultTranspose_Shell"
246 PetscErrorCode MatMultTranspose_Shell(Mat A,Vec x,Vec y)
247 {
248   Mat_Shell      *shell = (Mat_Shell*)A->data;
249   PetscErrorCode ierr;
250   Vec            xx;
251 
252   PetscFunctionBegin;
253   ierr = MatShellPreScaleLeft(A,x,&xx);CHKERRQ(ierr);
254   ierr = (*shell->multtranspose)(A,xx,y);CHKERRQ(ierr);
255   ierr = MatShellShiftAndScale(A,xx,y);CHKERRQ(ierr);
256   ierr = MatShellPostScaleRight(A,y);CHKERRQ(ierr);
257   PetscFunctionReturn(0);
258 }
259 
260 #undef __FUNCT__
261 #define __FUNCT__ "MatGetDiagonal_Shell"
262 PetscErrorCode MatGetDiagonal_Shell(Mat A,Vec v)
263 {
264   Mat_Shell      *shell = (Mat_Shell*)A->data;
265   PetscErrorCode ierr;
266 
267   PetscFunctionBegin;
268   ierr = (*shell->getdiagonal)(A,v);CHKERRQ(ierr);
269   ierr = VecScale(v,shell->vscale);CHKERRQ(ierr);
270   if (shell->dshift) {
271     ierr = VecPointwiseMult(v,v,shell->dshift);CHKERRQ(ierr);
272   } else {
273     ierr = VecShift(v,shell->vshift);CHKERRQ(ierr);
274   }
275   if (shell->left)  {ierr = VecPointwiseMult(v,v,shell->left);CHKERRQ(ierr);}
276   if (shell->right) {ierr = VecPointwiseMult(v,v,shell->right);CHKERRQ(ierr);}
277   PetscFunctionReturn(0);
278 }
279 
280 #undef __FUNCT__
281 #define __FUNCT__ "MatShift_Shell"
282 PetscErrorCode MatShift_Shell(Mat Y,PetscScalar a)
283 {
284   Mat_Shell *shell = (Mat_Shell*)Y->data;
285   PetscErrorCode ierr;
286 
287   PetscFunctionBegin;
288   if (shell->left || shell->right || shell->dshift) {
289     if (!shell->dshift) {
290       if (!shell->dshift_owned) {ierr = VecDuplicate(shell->left ? shell->left : shell->right, &shell->dshift_owned);CHKERRQ(ierr);}
291       shell->dshift = shell->dshift_owned;
292       ierr = VecSet(shell->dshift,shell->vshift+a);CHKERRQ(ierr);
293     } else {ierr = VecScale(shell->dshift,a);CHKERRQ(ierr);}
294     if (shell->left)  {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->left);CHKERRQ(ierr);}
295     if (shell->right) {ierr = VecPointwiseDivide(shell->dshift,shell->dshift,shell->right);CHKERRQ(ierr);}
296   } else shell->vshift += a;
297   ierr = MatShellUseScaledMethods(Y);CHKERRQ(ierr);
298   PetscFunctionReturn(0);
299 }
300 
301 #undef __FUNCT__
302 #define __FUNCT__ "MatScale_Shell"
303 PetscErrorCode MatScale_Shell(Mat Y,PetscScalar a)
304 {
305   Mat_Shell      *shell = (Mat_Shell*)Y->data;
306   PetscErrorCode ierr;
307 
308   PetscFunctionBegin;
309   shell->vscale *= a;
310   if (shell->dshift) {ierr = VecScale(shell->dshift,a);CHKERRQ(ierr);}
311   else shell->vshift *= a;
312   ierr = MatShellUseScaledMethods(Y);CHKERRQ(ierr);
313   PetscFunctionReturn(0);
314 }
315 
316 #undef __FUNCT__
317 #define __FUNCT__ "MatDiagonalScale_Shell"
318 static PetscErrorCode MatDiagonalScale_Shell(Mat Y,Vec left,Vec right)
319 {
320   Mat_Shell *shell = (Mat_Shell*)Y->data;
321   PetscErrorCode ierr;
322 
323   PetscFunctionBegin;
324   if (left) {
325     if (!shell->left_owned) {ierr = VecDuplicate(left,&shell->left_owned);CHKERRQ(ierr);}
326     if (shell->left) {ierr = VecPointwiseMult(shell->left,shell->left,left);CHKERRQ(ierr);}
327     else {
328       shell->left = shell->left_owned;
329       ierr = VecCopy(left,shell->left);CHKERRQ(ierr);
330     }
331   }
332   if (right) {
333     if (!shell->right_owned) {ierr = VecDuplicate(right,&shell->right_owned);CHKERRQ(ierr);}
334     if (shell->right) {ierr = VecPointwiseMult(shell->right,shell->right,right);CHKERRQ(ierr);}
335     else {
336       shell->right = shell->right_owned;
337       ierr = VecCopy(right,shell->right);CHKERRQ(ierr);
338     }
339   }
340   ierr = MatShellUseScaledMethods(Y);CHKERRQ(ierr);
341   PetscFunctionReturn(0);
342 }
343 
344 #undef __FUNCT__
345 #define __FUNCT__ "MatAssemblyEnd_Shell"
346 PetscErrorCode MatAssemblyEnd_Shell(Mat Y,MatAssemblyType t)
347 {
348   Mat_Shell *shell = (Mat_Shell*)Y->data;
349 
350   PetscFunctionBegin;
351   if (t == MAT_FINAL_ASSEMBLY) {
352     shell->vshift      = 0.0;
353     shell->vscale      = 1.0;
354     shell->dshift      = PETSC_NULL;
355     shell->left        = PETSC_NULL;
356     shell->right       = PETSC_NULL;
357     if (shell->mult) {
358       Y->ops->mult = shell->mult;
359       shell->mult  = PETSC_NULL;
360     }
361     if (shell->multtranspose) {
362       Y->ops->multtranspose = shell->multtranspose;
363       shell->multtranspose  = PETSC_NULL;
364     }
365     if (shell->getdiagonal) {
366       Y->ops->getdiagonal = shell->getdiagonal;
367       shell->getdiagonal  = PETSC_NULL;
368     }
369     shell->usingscaled = PETSC_FALSE;
370   }
371   PetscFunctionReturn(0);
372 }
373 
374 extern PetscErrorCode MatConvert_Shell(Mat, const MatType,MatReuse,Mat*);
375 
376 static struct _MatOps MatOps_Values = {0,
377        0,
378        0,
379        0,
380 /* 4*/ 0,
381        0,
382        0,
383        0,
384        0,
385        0,
386 /*10*/ 0,
387        0,
388        0,
389        0,
390        0,
391 /*15*/ 0,
392        0,
393        0,
394        MatDiagonalScale_Shell,
395        0,
396 /*20*/ 0,
397        MatAssemblyEnd_Shell,
398        0,
399        0,
400 /*24*/ 0,
401        0,
402        0,
403        0,
404        0,
405 /*29*/ 0,
406        0,
407        0,
408        0,
409        0,
410 /*34*/ 0,
411        0,
412        0,
413        0,
414        0,
415 /*39*/ 0,
416        0,
417        0,
418        0,
419        0,
420 /*44*/ 0,
421        MatScale_Shell,
422        MatShift_Shell,
423        0,
424        0,
425 /*49*/ 0,
426        0,
427        0,
428        0,
429        0,
430 /*54*/ 0,
431        0,
432        0,
433        0,
434        0,
435 /*59*/ 0,
436        MatDestroy_Shell,
437        0,
438        0,
439        0,
440 /*64*/ 0,
441        0,
442        0,
443        0,
444        0,
445 /*69*/ 0,
446        0,
447        MatConvert_Shell,
448        0,
449        0,
450 /*74*/ 0,
451        0,
452        0,
453        0,
454        0,
455 /*79*/ 0,
456        0,
457        0,
458        0,
459        0,
460 /*84*/ 0,
461        0,
462        0,
463        0,
464        0,
465 /*89*/ 0,
466        0,
467        0,
468        0,
469        0,
470 /*94*/ 0,
471        0,
472        0,
473        0};
474 
475 /*MC
476    MATSHELL - MATSHELL = "shell" - A matrix type to be used to define your own matrix type -- perhaps matrix free.
477 
478   Level: advanced
479 
480 .seealso: MatCreateShell
481 M*/
482 
483 EXTERN_C_BEGIN
484 #undef __FUNCT__
485 #define __FUNCT__ "MatCreate_Shell"
486 PetscErrorCode  MatCreate_Shell(Mat A)
487 {
488   Mat_Shell      *b;
489   PetscErrorCode ierr;
490 
491   PetscFunctionBegin;
492   ierr = PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
493 
494   ierr = PetscNewLog(A,Mat_Shell,&b);CHKERRQ(ierr);
495   A->data = (void*)b;
496 
497   ierr = PetscLayoutSetBlockSize(A->rmap,1);CHKERRQ(ierr);
498   ierr = PetscLayoutSetBlockSize(A->cmap,1);CHKERRQ(ierr);
499   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
500   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
501 
502   b->ctx           = 0;
503   b->vshift        = 0.0;
504   b->vscale        = 1.0;
505   b->mult          = 0;
506   b->multtranspose = 0;
507   b->getdiagonal   = 0;
508   A->assembled     = PETSC_TRUE;
509   A->preallocated  = PETSC_FALSE;
510   ierr = PetscObjectChangeTypeName((PetscObject)A,MATSHELL);CHKERRQ(ierr);
511   PetscFunctionReturn(0);
512 }
513 EXTERN_C_END
514 
515 #undef __FUNCT__
516 #define __FUNCT__ "MatCreateShell"
517 /*@C
518    MatCreateShell - Creates a new matrix class for use with a user-defined
519    private data storage format.
520 
521   Collective on MPI_Comm
522 
523    Input Parameters:
524 +  comm - MPI communicator
525 .  m - number of local rows (must be given)
526 .  n - number of local columns (must be given)
527 .  M - number of global rows (may be PETSC_DETERMINE)
528 .  N - number of global columns (may be PETSC_DETERMINE)
529 -  ctx - pointer to data needed by the shell matrix routines
530 
531    Output Parameter:
532 .  A - the matrix
533 
534    Level: advanced
535 
536   Usage:
537 $    extern int mult(Mat,Vec,Vec);
538 $    MatCreateShell(comm,m,n,M,N,ctx,&mat);
539 $    MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult);
540 $    [ Use matrix for operations that have been set ]
541 $    MatDestroy(mat);
542 
543    Notes:
544    The shell matrix type is intended to provide a simple class to use
545    with KSP (such as, for use with matrix-free methods). You should not
546    use the shell type if you plan to define a complete matrix class.
547 
548    Fortran Notes: The context can only be an integer or a PetscObject
549       unfortunately it cannot be a Fortran array or derived type.
550 
551    PETSc requires that matrices and vectors being used for certain
552    operations are partitioned accordingly.  For example, when
553    creating a shell matrix, A, that supports parallel matrix-vector
554    products using MatMult(A,x,y) the user should set the number
555    of local matrix rows to be the number of local elements of the
556    corresponding result vector, y. Note that this is information is
557    required for use of the matrix interface routines, even though
558    the shell matrix may not actually be physically partitioned.
559    For example,
560 
561 $
562 $     Vec x, y
563 $     extern int mult(Mat,Vec,Vec);
564 $     Mat A
565 $
566 $     VecCreateMPI(comm,PETSC_DECIDE,M,&y);
567 $     VecCreateMPI(comm,PETSC_DECIDE,N,&x);
568 $     VecGetLocalSize(y,&m);
569 $     VecGetLocalSize(x,&n);
570 $     MatCreateShell(comm,m,n,M,N,ctx,&A);
571 $     MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult);
572 $     MatMult(A,x,y);
573 $     MatDestroy(A);
574 $     VecDestroy(y); VecDestroy(x);
575 $
576 
577 .keywords: matrix, shell, create
578 
579 .seealso: MatShellSetOperation(), MatHasOperation(), MatShellGetContext(), MatShellSetContext()
580 @*/
581 PetscErrorCode  MatCreateShell(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx,Mat *A)
582 {
583   PetscErrorCode ierr;
584 
585   PetscFunctionBegin;
586   ierr = MatCreate(comm,A);CHKERRQ(ierr);
587   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
588   ierr = MatSetType(*A,MATSHELL);CHKERRQ(ierr);
589   ierr = MatShellSetContext(*A,ctx);CHKERRQ(ierr);
590   ierr = MatSetUp(*A);CHKERRQ(ierr);
591   PetscFunctionReturn(0);
592 }
593 
594 #undef __FUNCT__
595 #define __FUNCT__ "MatShellSetContext"
596 /*@
597     MatShellSetContext - sets the context for a shell matrix
598 
599    Logically Collective on Mat
600 
601     Input Parameters:
602 +   mat - the shell matrix
603 -   ctx - the context
604 
605    Level: advanced
606 
607    Fortran Notes: The context can only be an integer or a PetscObject
608       unfortunately it cannot be a Fortran array or derived type.
609 
610 .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation()
611 @*/
612 PetscErrorCode  MatShellSetContext(Mat mat,void *ctx)
613 {
614   Mat_Shell      *shell = (Mat_Shell*)mat->data;
615   PetscErrorCode ierr;
616   PetscBool      flg;
617 
618   PetscFunctionBegin;
619   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
620   ierr = PetscTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
621   if (flg) {
622     shell->ctx = ctx;
623   } else SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_SUP,"Cannot attach context to non-shell matrix");
624   PetscFunctionReturn(0);
625 }
626 
627 #undef __FUNCT__
628 #define __FUNCT__ "MatShellSetOperation"
629 /*@C
630     MatShellSetOperation - Allows user to set a matrix operation for
631                            a shell matrix.
632 
633    Logically Collective on Mat
634 
635     Input Parameters:
636 +   mat - the shell matrix
637 .   op - the name of the operation
638 -   f - the function that provides the operation.
639 
640    Level: advanced
641 
642     Usage:
643 $      extern PetscErrorCode usermult(Mat,Vec,Vec);
644 $      ierr = MatCreateShell(comm,m,n,M,N,ctx,&A);
645 $      ierr = MatShellSetOperation(A,MATOP_MULT,(void(*)(void))usermult);
646 
647     Notes:
648     See the file include/petscmat.h for a complete list of matrix
649     operations, which all have the form MATOP_<OPERATION>, where
650     <OPERATION> is the name (in all capital letters) of the
651     user interface routine (e.g., MatMult() -> MATOP_MULT).
652 
653     All user-provided functions (execept for MATOP_DESTROY) should have the same calling
654     sequence as the usual matrix interface routines, since they
655     are intended to be accessed via the usual matrix interface
656     routines, e.g.,
657 $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)
658 
659     In particular each function MUST return an error code of 0 on success and
660     nonzero on failure.
661 
662     Within each user-defined routine, the user should call
663     MatShellGetContext() to obtain the user-defined context that was
664     set by MatCreateShell().
665 
666     Fortran Notes: For MatGetVecs() the user code should check if the input left or right matrix is -1 and in that case not
667        generate a matrix. See src/mat/examples/tests/ex120f.F
668 
669 .keywords: matrix, shell, set, operation
670 
671 .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext()
672 @*/
673 PetscErrorCode  MatShellSetOperation(Mat mat,MatOperation op,void (*f)(void))
674 {
675   PetscErrorCode ierr;
676   PetscBool      flg;
677 
678   PetscFunctionBegin;
679   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
680   if (op == MATOP_DESTROY) {
681     ierr = PetscTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
682     if (flg) {
683        Mat_Shell *shell = (Mat_Shell*)mat->data;
684        shell->destroy                 = (PetscErrorCode (*)(Mat)) f;
685     } else mat->ops->destroy          = (PetscErrorCode (*)(Mat)) f;
686   }
687   else if (op == MATOP_VIEW) mat->ops->view  = (PetscErrorCode (*)(Mat,PetscViewer)) f;
688   else                       (((void(**)(void))mat->ops)[op]) = f;
689 
690   PetscFunctionReturn(0);
691 }
692 
693 #undef __FUNCT__
694 #define __FUNCT__ "MatShellGetOperation"
695 /*@C
696     MatShellGetOperation - Gets a matrix function for a shell matrix.
697 
698     Not Collective
699 
700     Input Parameters:
701 +   mat - the shell matrix
702 -   op - the name of the operation
703 
704     Output Parameter:
705 .   f - the function that provides the operation.
706 
707     Level: advanced
708 
709     Notes:
710     See the file include/petscmat.h for a complete list of matrix
711     operations, which all have the form MATOP_<OPERATION>, where
712     <OPERATION> is the name (in all capital letters) of the
713     user interface routine (e.g., MatMult() -> MATOP_MULT).
714 
715     All user-provided functions have the same calling
716     sequence as the usual matrix interface routines, since they
717     are intended to be accessed via the usual matrix interface
718     routines, e.g.,
719 $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)
720 
721     Within each user-defined routine, the user should call
722     MatShellGetContext() to obtain the user-defined context that was
723     set by MatCreateShell().
724 
725 .keywords: matrix, shell, set, operation
726 
727 .seealso: MatCreateShell(), MatShellGetContext(), MatShellSetOperation(), MatShellSetContext()
728 @*/
729 PetscErrorCode  MatShellGetOperation(Mat mat,MatOperation op,void(**f)(void))
730 {
731   PetscErrorCode ierr;
732   PetscBool      flg;
733 
734   PetscFunctionBegin;
735   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
736   if (op == MATOP_DESTROY) {
737     ierr = PetscTypeCompare((PetscObject)mat,MATSHELL,&flg);CHKERRQ(ierr);
738     if (flg) {
739       Mat_Shell *shell = (Mat_Shell*)mat->data;
740       *f = (void(*)(void))shell->destroy;
741     } else {
742       *f = (void(*)(void))mat->ops->destroy;
743     }
744   } else if (op == MATOP_VIEW) {
745     *f = (void(*)(void))mat->ops->view;
746   } else {
747     *f = (((void(**)(void))mat->ops)[op]);
748   }
749 
750   PetscFunctionReturn(0);
751 }
752 
753