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