Actual source code: shell.c

petsc-master 2019-05-21
Report Typos and Errors

  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>
  9:  #include <petsc/private/vecimpl.h>

 11: struct _MatShellOps {
 12:   /*  3 */ PetscErrorCode (*mult)(Mat,Vec,Vec);
 13:   /*  5 */ PetscErrorCode (*multtranspose)(Mat,Vec,Vec);
 14:   /* 17 */ PetscErrorCode (*getdiagonal)(Mat,Vec);
 15:   /* 43 */ PetscErrorCode (*copy)(Mat,Mat,MatStructure);
 16:   /* 60 */ PetscErrorCode (*destroy)(Mat);
 17: };

 19: typedef struct {
 20:   struct _MatShellOps ops[1];

 22:   PetscScalar vscale,vshift;
 23:   Vec         dshift;
 24:   Vec         left,right;
 25:   Vec         left_work,right_work;
 26:   Vec         left_add_work,right_add_work;
 27:   Mat         axpy;
 28:   PetscScalar axpy_vscale;
 29:   PetscBool   managescalingshifts;                   /* The user will manage the scaling and shifts for the MATSHELL, not the default */
 30:   /* support for ZeroRows/Columns operations */
 31:   IS          zrows;
 32:   IS          zcols;
 33:   Vec         zvals;
 34:   Vec         zvals_w;
 35:   VecScatter  zvals_sct_r;
 36:   VecScatter  zvals_sct_c;
 37:   void        *ctx;
 38: } Mat_Shell;


 41: /*
 42:      Store and scale values on zeroed rows
 43:      xx = [x_1, 0], 0 on zeroed columns
 44: */
 45: static PetscErrorCode MatShellPreZeroRight(Mat A,Vec x,Vec *xx)
 46: {
 47:   Mat_Shell      *shell = (Mat_Shell*)A->data;

 51:   *xx = x;
 52:   if (shell->zrows) {
 53:     VecSet(shell->zvals_w,0.0);
 54:     VecScatterBegin(shell->zvals_sct_c,x,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD);
 55:     VecScatterEnd(shell->zvals_sct_c,x,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD);
 56:     VecPointwiseMult(shell->zvals_w,shell->zvals_w,shell->zvals);
 57:   }
 58:   if (shell->zcols) {
 59:     if (!shell->right_work) {
 60:       MatCreateVecs(A,&shell->right_work,NULL);
 61:     }
 62:     VecCopy(x,shell->right_work);
 63:     VecISSet(shell->right_work,shell->zcols,0.0);
 64:     *xx  = shell->right_work;
 65:   }
 66:   return(0);
 67: }

 69: /* Insert properly diagonally scaled values stored in MatShellPreZeroRight */
 70: static PetscErrorCode MatShellPostZeroLeft(Mat A,Vec x)
 71: {
 72:   Mat_Shell      *shell = (Mat_Shell*)A->data;

 76:   if (shell->zrows) {
 77:     VecScatterBegin(shell->zvals_sct_r,shell->zvals_w,x,INSERT_VALUES,SCATTER_REVERSE);
 78:     VecScatterEnd(shell->zvals_sct_r,shell->zvals_w,x,INSERT_VALUES,SCATTER_REVERSE);
 79:   }
 80:   return(0);
 81: }

 83: /*
 84:      Store and scale values on zeroed rows
 85:      xx = [x_1, 0], 0 on zeroed rows
 86: */
 87: static PetscErrorCode MatShellPreZeroLeft(Mat A,Vec x,Vec *xx)
 88: {
 89:   Mat_Shell      *shell = (Mat_Shell*)A->data;

 93:   *xx = NULL;
 94:   if (!shell->zrows) {
 95:     *xx = x;
 96:   } else {
 97:     if (!shell->left_work) {
 98:       MatCreateVecs(A,NULL,&shell->left_work);
 99:     }
100:     VecCopy(x,shell->left_work);
101:     VecSet(shell->zvals_w,0.0);
102:     VecScatterBegin(shell->zvals_sct_r,shell->zvals_w,shell->left_work,INSERT_VALUES,SCATTER_REVERSE);
103:     VecScatterEnd(shell->zvals_sct_r,shell->zvals_w,shell->left_work,INSERT_VALUES,SCATTER_REVERSE);
104:     VecScatterBegin(shell->zvals_sct_r,x,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD);
105:     VecScatterEnd(shell->zvals_sct_r,x,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD);
106:     VecPointwiseMult(shell->zvals_w,shell->zvals_w,shell->zvals);
107:     *xx  = shell->left_work;
108:   }
109:   return(0);
110: }

112: /* Zero zero-columns contributions, sum contributions from properly scaled values stored in MatShellPreZeroLeft */
113: static PetscErrorCode MatShellPostZeroRight(Mat A,Vec x)
114: {
115:   Mat_Shell      *shell = (Mat_Shell*)A->data;

119:   if (shell->zcols) {
120:     VecISSet(x,shell->zcols,0.0);
121:   }
122:   if (shell->zrows) {
123:     VecScatterBegin(shell->zvals_sct_c,shell->zvals_w,x,ADD_VALUES,SCATTER_REVERSE);
124:     VecScatterEnd(shell->zvals_sct_c,shell->zvals_w,x,ADD_VALUES,SCATTER_REVERSE);
125:   }
126:   return(0);
127: }

129: /*
130:       xx = diag(left)*x
131: */
132: static PetscErrorCode MatShellPreScaleLeft(Mat A,Vec x,Vec *xx)
133: {
134:   Mat_Shell      *shell = (Mat_Shell*)A->data;

138:   *xx = NULL;
139:   if (!shell->left) {
140:     *xx = x;
141:   } else {
142:     if (!shell->left_work) {VecDuplicate(shell->left,&shell->left_work);}
143:     VecPointwiseMult(shell->left_work,x,shell->left);
144:     *xx  = shell->left_work;
145:   }
146:   return(0);
147: }

149: /*
150:      xx = diag(right)*x
151: */
152: static PetscErrorCode MatShellPreScaleRight(Mat A,Vec x,Vec *xx)
153: {
154:   Mat_Shell      *shell = (Mat_Shell*)A->data;

158:   *xx = NULL;
159:   if (!shell->right) {
160:     *xx = x;
161:   } else {
162:     if (!shell->right_work) {VecDuplicate(shell->right,&shell->right_work);}
163:     VecPointwiseMult(shell->right_work,x,shell->right);
164:     *xx  = shell->right_work;
165:   }
166:   return(0);
167: }

169: /*
170:     x = diag(left)*x
171: */
172: static PetscErrorCode MatShellPostScaleLeft(Mat A,Vec x)
173: {
174:   Mat_Shell      *shell = (Mat_Shell*)A->data;

178:   if (shell->left) {VecPointwiseMult(x,x,shell->left);}
179:   return(0);
180: }

182: /*
183:     x = diag(right)*x
184: */
185: static PetscErrorCode MatShellPostScaleRight(Mat A,Vec x)
186: {
187:   Mat_Shell      *shell = (Mat_Shell*)A->data;

191:   if (shell->right) {VecPointwiseMult(x,x,shell->right);}
192:   return(0);
193: }

195: /*
196:          Y = vscale*Y + diag(dshift)*X + vshift*X

198:          On input Y already contains A*x
199: */
200: static PetscErrorCode MatShellShiftAndScale(Mat A,Vec X,Vec Y)
201: {
202:   Mat_Shell      *shell = (Mat_Shell*)A->data;

206:   if (shell->dshift) {          /* get arrays because there is no VecPointwiseMultAdd() */
207:     PetscInt          i,m;
208:     const PetscScalar *x,*d;
209:     PetscScalar       *y;
210:     VecGetLocalSize(X,&m);
211:     VecGetArrayRead(shell->dshift,&d);
212:     VecGetArrayRead(X,&x);
213:     VecGetArray(Y,&y);
214:     for (i=0; i<m; i++) y[i] = shell->vscale*y[i] + d[i]*x[i];
215:     VecRestoreArrayRead(shell->dshift,&d);
216:     VecRestoreArrayRead(X,&x);
217:     VecRestoreArray(Y,&y);
218:   } else {
219:     VecScale(Y,shell->vscale);
220:   }
221:   if (shell->vshift != 0.0) {VecAXPY(Y,shell->vshift,X);} /* if test is for non-square matrices */
222:   return(0);
223: }

225: /*@
226:     MatShellGetContext - Returns the user-provided context associated with a shell matrix.

228:     Not Collective

230:     Input Parameter:
231: .   mat - the matrix, should have been created with MatCreateShell()

233:     Output Parameter:
234: .   ctx - the user provided context

236:     Level: advanced

238:    Fortran Notes:
239:     To use this from Fortran you must write a Fortran interface definition for this
240:     function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.

242: .keywords: matrix, shell, get, context

244: .seealso: MatCreateShell(), MatShellSetOperation(), MatShellSetContext()
245: @*/
246: PetscErrorCode  MatShellGetContext(Mat mat,void *ctx)
247: {
249:   PetscBool      flg;

254:   PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);
255:   if (flg) *(void**)ctx = ((Mat_Shell*)(mat->data))->ctx;
256:   else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot get context from non-shell matrix");
257:   return(0);
258: }

260: static PetscErrorCode MatZeroRowsColumns_Local_Shell(Mat mat,PetscInt nr,PetscInt rows[],PetscInt nc,PetscInt cols[],PetscScalar diag,PetscBool rc)
261: {
263:   Mat_Shell      *shell = (Mat_Shell*)mat->data;
264:   Vec            x = NULL,b = NULL;
265:   IS             is1, is2;
266:   const PetscInt *ridxs;
267:   PetscInt       *idxs,*gidxs;
268:   PetscInt       cum,rst,cst,i;

271:   if (shell->axpy) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot diagonal scale MATSHELL after MatAXPY operation");
272:   if (!shell->zvals) {
273:     MatCreateVecs(mat,NULL,&shell->zvals);
274:   }
275:   if (!shell->zvals_w) {
276:     VecDuplicate(shell->zvals,&shell->zvals_w);
277:   }
278:   MatGetOwnershipRange(mat,&rst,NULL);
279:   MatGetOwnershipRangeColumn(mat,&cst,NULL);

281:   /* Expand/create index set of zeroed rows */
282:   PetscMalloc1(nr,&idxs);
283:   for (i = 0; i < nr; i++) idxs[i] = rows[i] + rst;
284:   ISCreateGeneral(PETSC_COMM_SELF,nr,idxs,PETSC_OWN_POINTER,&is1);
285:   ISSort(is1);
286:   VecISSet(shell->zvals,is1,diag);
287:   if (shell->zrows) {
288:     ISSum(shell->zrows,is1,&is2);
289:     ISDestroy(&shell->zrows);
290:     ISDestroy(&is1);
291:     shell->zrows = is2;
292:   } else shell->zrows = is1;

294:   /* Create scatters for diagonal values communications */
295:   VecScatterDestroy(&shell->zvals_sct_c);
296:   VecScatterDestroy(&shell->zvals_sct_r);

298:   /* row scatter: from/to left vector */
299:   MatCreateVecs(mat,&x,&b);
300:   VecScatterCreate(b,shell->zrows,shell->zvals_w,shell->zrows,&shell->zvals_sct_r);

302:   /* col scatter: from right vector to left vector */
303:   ISGetIndices(shell->zrows,&ridxs);
304:   ISGetLocalSize(shell->zrows,&nr);
305:   PetscMalloc1(nr,&gidxs);
306:   for (i = 0, cum  = 0; i < nr; i++) {
307:     if (ridxs[i] >= mat->cmap->N) continue;
308:     gidxs[cum] = ridxs[i];
309:     cum++;
310:   }
311:   ISCreateGeneral(PetscObjectComm((PetscObject)mat),cum,gidxs,PETSC_OWN_POINTER,&is1);
312:   VecScatterCreate(x,is1,shell->zvals_w,is1,&shell->zvals_sct_c);
313:   ISDestroy(&is1);
314:   VecDestroy(&x);
315:   VecDestroy(&b);

317:   /* Expand/create index set of zeroed columns */
318:   if (rc) {
319:     PetscMalloc1(nc,&idxs);
320:     for (i = 0; i < nc; i++) idxs[i] = cols[i] + cst;
321:     ISCreateGeneral(PETSC_COMM_SELF,nc,idxs,PETSC_OWN_POINTER,&is1);
322:     ISSort(is1);
323:     if (shell->zcols) {
324:       ISSum(shell->zcols,is1,&is2);
325:       ISDestroy(&shell->zcols);
326:       ISDestroy(&is1);
327:       shell->zcols = is2;
328:     } else shell->zcols = is1;
329:   }
330:   return(0);
331: }

333: static PetscErrorCode MatZeroRows_Shell(Mat mat,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
334: {
335:   PetscInt       nr, *lrows;

339:   if (x && b) {
340:     Vec          xt;
341:     PetscScalar *vals;
342:     PetscInt    *gcols,i,st,nl,nc;

344:     PetscMalloc1(n,&gcols);
345:     for (i = 0, nc = 0; i < n; i++) if (rows[i] < mat->cmap->N) gcols[nc++] = rows[i];

347:     MatCreateVecs(mat,&xt,NULL);
348:     VecCopy(x,xt);
349:     PetscCalloc1(nc,&vals);
350:     VecSetValues(xt,nc,gcols,vals,INSERT_VALUES); /* xt = [x1, 0] */
351:     PetscFree(vals);
352:     VecAssemblyBegin(xt);
353:     VecAssemblyEnd(xt);
354:     VecAYPX(xt,-1.0,x);                           /* xt = [0, x2] */

356:     VecGetOwnershipRange(xt,&st,NULL);
357:     VecGetLocalSize(xt,&nl);
358:     VecGetArray(xt,&vals);
359:     for (i = 0; i < nl; i++) {
360:       PetscInt g = i + st;
361:       if (g > mat->rmap->N) continue;
362:       if (PetscAbsScalar(vals[i]) == 0.0) continue;
363:       VecSetValue(b,g,diag*vals[i],INSERT_VALUES);
364:     }
365:     VecRestoreArray(xt,&vals);
366:     VecAssemblyBegin(b);
367:     VecAssemblyEnd(b);                            /* b  = [b1, x2 * diag] */
368:     VecDestroy(&xt);
369:     PetscFree(gcols);
370:   }
371:   PetscLayoutMapLocal(mat->rmap,n,rows,&nr,&lrows,NULL);
372:   MatZeroRowsColumns_Local_Shell(mat,nr,lrows,0,NULL,diag,PETSC_FALSE);
373:   PetscFree(lrows);
374:   return(0);
375: }

377: static PetscErrorCode MatZeroRowsColumns_Shell(Mat mat,PetscInt n,const PetscInt rowscols[],PetscScalar diag,Vec x,Vec b)
378: {
379:   PetscInt       *lrows, *lcols;
380:   PetscInt       nr, nc;
381:   PetscBool      congruent;

385:   if (x && b) {
386:     Vec          xt, bt;
387:     PetscScalar *vals;
388:     PetscInt    *grows,*gcols,i,st,nl;

390:     PetscMalloc2(n,&grows,n,&gcols);
391:     for (i = 0, nr = 0; i < n; i++) if (rowscols[i] < mat->rmap->N) grows[nr++] = rowscols[i];
392:     for (i = 0, nc = 0; i < n; i++) if (rowscols[i] < mat->cmap->N) gcols[nc++] = rowscols[i];
393:     PetscCalloc1(n,&vals);

395:     MatCreateVecs(mat,&xt,&bt);
396:     VecCopy(x,xt);
397:     VecSetValues(xt,nc,gcols,vals,INSERT_VALUES); /* xt = [x1, 0] */
398:     VecAssemblyBegin(xt);
399:     VecAssemblyEnd(xt);
400:     VecAXPY(xt,-1.0,x);                           /* xt = [0, -x2] */
401:     MatMult(mat,xt,bt);                           /* bt = [-A12*x2,-A22*x2] */
402:     VecSetValues(bt,nr,grows,vals,INSERT_VALUES); /* bt = [-A12*x2,0] */
403:     VecAssemblyBegin(bt);
404:     VecAssemblyEnd(bt);
405:     VecAXPY(b,1.0,bt);                            /* b  = [b1 - A12*x2, b2] */
406:     VecSetValues(bt,nr,grows,vals,INSERT_VALUES); /* b  = [b1 - A12*x2, 0] */
407:     VecAssemblyBegin(bt);
408:     VecAssemblyEnd(bt);
409:     PetscFree(vals);

411:     VecGetOwnershipRange(xt,&st,NULL);
412:     VecGetLocalSize(xt,&nl);
413:     VecGetArray(xt,&vals);
414:     for (i = 0; i < nl; i++) {
415:       PetscInt g = i + st;
416:       if (g > mat->rmap->N) continue;
417:       if (PetscAbsScalar(vals[i]) == 0.0) continue;
418:       VecSetValue(b,g,-diag*vals[i],INSERT_VALUES);
419:     }
420:     VecRestoreArray(xt,&vals);
421:     VecAssemblyBegin(b);
422:     VecAssemblyEnd(b);                            /* b  = [b1 - A12*x2, x2 * diag] */
423:     VecDestroy(&xt);
424:     VecDestroy(&bt);
425:     PetscFree2(grows,gcols);
426:   }
427:   PetscLayoutMapLocal(mat->rmap,n,rowscols,&nr,&lrows,NULL);
428:   MatHasCongruentLayouts(mat,&congruent);
429:   if (congruent) {
430:     nc    = nr;
431:     lcols = lrows;
432:   } else { /* MatZeroRowsColumns implicitly assumes the rowscols indices are for a square matrix, here we handle a more general case */
433:     PetscInt i,nt,*t;

435:     PetscMalloc1(n,&t);
436:     for (i = 0, nt = 0; i < n; i++) if (rowscols[i] < mat->cmap->N) t[nt++] = rowscols[i];
437:     PetscLayoutMapLocal(mat->cmap,nt,t,&nc,&lcols,NULL);
438:     PetscFree(t);
439:   }
440:   MatZeroRowsColumns_Local_Shell(mat,nr,lrows,nc,lcols,diag,PETSC_TRUE);
441:   if (!congruent) {
442:     PetscFree(lcols);
443:   }
444:   PetscFree(lrows);
445:   return(0);
446: }

448: PetscErrorCode MatDestroy_Shell(Mat mat)
449: {
451:   Mat_Shell      *shell = (Mat_Shell*)mat->data;

454:   if (shell->ops->destroy) {
455:     (*shell->ops->destroy)(mat);
456:   }
457:   VecDestroy(&shell->left);
458:   VecDestroy(&shell->right);
459:   VecDestroy(&shell->dshift);
460:   VecDestroy(&shell->left_work);
461:   VecDestroy(&shell->right_work);
462:   VecDestroy(&shell->left_add_work);
463:   VecDestroy(&shell->right_add_work);
464:   MatDestroy(&shell->axpy);
465:   VecDestroy(&shell->zvals_w);
466:   VecDestroy(&shell->zvals);
467:   VecScatterDestroy(&shell->zvals_sct_c);
468:   VecScatterDestroy(&shell->zvals_sct_r);
469:   ISDestroy(&shell->zrows);
470:   ISDestroy(&shell->zcols);
471:   PetscFree(mat->data);
472:   return(0);
473: }

475: PetscErrorCode MatCopy_Shell(Mat A,Mat B,MatStructure str)
476: {
477:   Mat_Shell       *shellA = (Mat_Shell*)A->data,*shellB = (Mat_Shell*)B->data;
478:   PetscErrorCode  ierr;
479:   PetscBool       matflg;

482:   PetscObjectTypeCompare((PetscObject)B,MATSHELL,&matflg);
483:   if (!matflg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_NOTSAMETYPE,"Output matrix must be a MATSHELL");

485:   PetscMemcpy(B->ops,A->ops,sizeof(struct _MatOps));
486:   PetscMemcpy(shellB->ops,shellA->ops,sizeof(struct _MatShellOps));

488:   if (shellA->ops->copy) {
489:     (*shellA->ops->copy)(A,B,str);
490:   }
491:   shellB->vscale = shellA->vscale;
492:   shellB->vshift = shellA->vshift;
493:   if (shellA->dshift) {
494:     if (!shellB->dshift) {
495:       VecDuplicate(shellA->dshift,&shellB->dshift);
496:     }
497:     VecCopy(shellA->dshift,shellB->dshift);
498:   } else {
499:     VecDestroy(&shellB->dshift);
500:   }
501:   if (shellA->left) {
502:     if (!shellB->left) {
503:       VecDuplicate(shellA->left,&shellB->left);
504:     }
505:     VecCopy(shellA->left,shellB->left);
506:   } else {
507:     VecDestroy(&shellB->left);
508:   }
509:   if (shellA->right) {
510:     if (!shellB->right) {
511:       VecDuplicate(shellA->right,&shellB->right);
512:     }
513:     VecCopy(shellA->right,shellB->right);
514:   } else {
515:     VecDestroy(&shellB->right);
516:   }
517:   MatDestroy(&shellB->axpy);
518:   if (shellA->axpy) {
519:     PetscObjectReference((PetscObject)shellA->axpy);
520:     shellB->axpy        = shellA->axpy;
521:     shellB->axpy_vscale = shellA->axpy_vscale;
522:   }
523:   if (shellA->zrows) {
524:     ISDuplicate(shellA->zrows,&shellB->zrows);
525:     if (shellA->zcols) {
526:       ISDuplicate(shellA->zcols,&shellB->zcols);
527:     }
528:     VecDuplicate(shellA->zvals,&shellB->zvals);
529:     VecCopy(shellA->zvals,shellB->zvals);
530:     VecDuplicate(shellA->zvals_w,&shellB->zvals_w);
531:     PetscObjectReference((PetscObject)shellA->zvals_sct_r);
532:     PetscObjectReference((PetscObject)shellA->zvals_sct_c);
533:     shellB->zvals_sct_r = shellA->zvals_sct_r;
534:     shellB->zvals_sct_c = shellA->zvals_sct_c;
535:   }
536:   return(0);
537: }

539: PetscErrorCode MatDuplicate_Shell(Mat mat,MatDuplicateOption op,Mat *M)
540: {
542:   void           *ctx;

545:   MatShellGetContext(mat,&ctx);
546:   MatCreateShell(PetscObjectComm((PetscObject)mat),mat->rmap->n,mat->cmap->n,mat->rmap->N,mat->cmap->N,ctx,M);
547:   if (op != MAT_DO_NOT_COPY_VALUES) {
548:     MatCopy(mat,*M,SAME_NONZERO_PATTERN);
549:   }
550:   return(0);
551: }

553: PetscErrorCode MatMult_Shell(Mat A,Vec x,Vec y)
554: {
555:   Mat_Shell        *shell = (Mat_Shell*)A->data;
556:   PetscErrorCode   ierr;
557:   Vec              xx;
558:   PetscObjectState instate,outstate;

561:   if (!shell->ops->mult) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Have not provided a MatMult() for this MATSHELL");
562:   MatShellPreZeroRight(A,x,&xx);
563:   MatShellPreScaleRight(A,xx,&xx);
564:   PetscObjectStateGet((PetscObject)y, &instate);
565:   (*shell->ops->mult)(A,xx,y);
566:   PetscObjectStateGet((PetscObject)y, &outstate);
567:   if (instate == outstate) {
568:     /* increase the state of the output vector since the user did not update its state themself as should have been done */
569:     PetscObjectStateIncrease((PetscObject)y);
570:   }
571:   MatShellShiftAndScale(A,xx,y);
572:   MatShellPostScaleLeft(A,y);
573:   MatShellPostZeroLeft(A,y);

575:   if (shell->axpy) {
576:     if (!shell->left_work) {MatCreateVecs(A,&shell->left_work,NULL);}
577:     MatMult(shell->axpy,x,shell->left_work);
578:     VecAXPY(y,shell->axpy_vscale,shell->left_work);
579:   }
580:   return(0);
581: }

583: PetscErrorCode MatMultAdd_Shell(Mat A,Vec x,Vec y,Vec z)
584: {
585:   Mat_Shell      *shell = (Mat_Shell*)A->data;

589:   if (y == z) {
590:     if (!shell->right_add_work) {VecDuplicate(z,&shell->right_add_work);}
591:     MatMult(A,x,shell->right_add_work);
592:     VecAXPY(z,1.0,shell->right_add_work);
593:   } else {
594:     MatMult(A,x,z);
595:     VecAXPY(z,1.0,y);
596:   }
597:   return(0);
598: }

600: PetscErrorCode MatMultTranspose_Shell(Mat A,Vec x,Vec y)
601: {
602:   Mat_Shell        *shell = (Mat_Shell*)A->data;
603:   PetscErrorCode   ierr;
604:   Vec              xx;
605:   PetscObjectState instate,outstate;

608:   MatShellPreZeroLeft(A,x,&xx);
609:   MatShellPreScaleLeft(A,xx,&xx);
610:   PetscObjectStateGet((PetscObject)y, &instate);
611:   if (!shell->ops->multtranspose) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Have not provided a MatMultTranspose() for this MATSHELL");
612:   (*shell->ops->multtranspose)(A,xx,y);
613:   PetscObjectStateGet((PetscObject)y, &outstate);
614:   if (instate == outstate) {
615:     /* increase the state of the output vector since the user did not update its state themself as should have been done */
616:     PetscObjectStateIncrease((PetscObject)y);
617:   }
618:   MatShellShiftAndScale(A,xx,y);
619:   MatShellPostScaleRight(A,y);
620:   MatShellPostZeroRight(A,y);

622:   if (shell->axpy) {
623:     if (!shell->right_work) {MatCreateVecs(A,NULL,&shell->right_work);}
624:     MatMultTranspose(shell->axpy,x,shell->right_work);
625:     VecAXPY(y,shell->axpy_vscale,shell->right_work);
626:   }
627:   return(0);
628: }

630: PetscErrorCode MatMultTransposeAdd_Shell(Mat A,Vec x,Vec y,Vec z)
631: {
632:   Mat_Shell      *shell = (Mat_Shell*)A->data;

636:   if (y == z) {
637:     if (!shell->left_add_work) {VecDuplicate(z,&shell->left_add_work);}
638:     MatMultTranspose(A,x,shell->left_add_work);
639:     VecAXPY(z,1.0,shell->left_add_work);
640:   } else {
641:     MatMultTranspose(A,x,z);
642:     VecAXPY(z,1.0,y);
643:   }
644:   return(0);
645: }

647: /*
648:           diag(left)(vscale*A + diag(dshift) + vshift I)diag(right)
649: */
650: PetscErrorCode MatGetDiagonal_Shell(Mat A,Vec v)
651: {
652:   Mat_Shell      *shell = (Mat_Shell*)A->data;

656:   if (shell->ops->getdiagonal) {
657:     (*shell->ops->getdiagonal)(A,v);
658:   } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Must provide shell matrix with routine to return diagonal using\nMatShellSetOperation(S,MATOP_GET_DIAGONAL,...)");
659:   VecScale(v,shell->vscale);
660:   if (shell->dshift) {
661:     VecAXPY(v,1.0,shell->dshift);
662:   }
663:   VecShift(v,shell->vshift);
664:   if (shell->left)  {VecPointwiseMult(v,v,shell->left);}
665:   if (shell->right) {VecPointwiseMult(v,v,shell->right);}
666:   if (shell->zrows) {
667:     VecScatterBegin(shell->zvals_sct_r,shell->zvals,v,INSERT_VALUES,SCATTER_REVERSE);
668:     VecScatterEnd(shell->zvals_sct_r,shell->zvals,v,INSERT_VALUES,SCATTER_REVERSE);
669:   }
670:   if (shell->axpy) {
671:     if (!shell->left_work) {VecDuplicate(v,&shell->left_work);}
672:     MatGetDiagonal(shell->axpy,shell->left_work);
673:     VecAXPY(v,shell->axpy_vscale,shell->left_work);
674:   }
675:   return(0);
676: }

678: PetscErrorCode MatShift_Shell(Mat Y,PetscScalar a)
679: {
680:   Mat_Shell      *shell = (Mat_Shell*)Y->data;

684:   if (shell->left || shell->right) {
685:     if (!shell->dshift) {
686:       VecDuplicate(shell->left ? shell->left : shell->right, &shell->dshift);
687:       VecSet(shell->dshift,a);
688:     } else {
689:       if (shell->left)  {VecPointwiseMult(shell->dshift,shell->dshift,shell->left);}
690:       if (shell->right) {VecPointwiseMult(shell->dshift,shell->dshift,shell->right);}
691:       VecShift(shell->dshift,a);
692:     }
693:     if (shell->left)  {VecPointwiseDivide(shell->dshift,shell->dshift,shell->left);}
694:     if (shell->right) {VecPointwiseDivide(shell->dshift,shell->dshift,shell->right);}
695:   } else shell->vshift += a;
696:   if (shell->zrows) {
697:     VecShift(shell->zvals,a);
698:   }
699:   return(0);
700: }

702: PetscErrorCode MatDiagonalSet_Shell_Private(Mat A,Vec D,PetscScalar s)
703: {
704:   Mat_Shell      *shell = (Mat_Shell*)A->data;

708:   if (!shell->dshift) {VecDuplicate(D,&shell->dshift);}
709:   if (shell->left || shell->right) {
710:     if (!shell->right_work) {VecDuplicate(shell->left ? shell->left : shell->right, &shell->right_work);}
711:     if (shell->left && shell->right)  {
712:       VecPointwiseDivide(shell->right_work,D,shell->left);
713:       VecPointwiseDivide(shell->right_work,shell->right_work,shell->right);
714:     } else if (shell->left) {
715:       VecPointwiseDivide(shell->right_work,D,shell->left);
716:     } else {
717:       VecPointwiseDivide(shell->right_work,D,shell->right);
718:     }
719:     VecAXPY(shell->dshift,s,shell->right_work);
720:   } else {
721:     VecAXPY(shell->dshift,s,D);
722:   }
723:   return(0);
724: }

726: PetscErrorCode MatDiagonalSet_Shell(Mat A,Vec D,InsertMode ins)
727: {
728:   Mat_Shell      *shell = (Mat_Shell*)A->data;
729:   Vec            d;

733:   if (ins == INSERT_VALUES) {
734:     if (!A->ops->getdiagonal) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Operation MATOP_GETDIAGONAL must be set first");
735:     VecDuplicate(D,&d);
736:     MatGetDiagonal(A,d);
737:     MatDiagonalSet_Shell_Private(A,d,-1.);
738:     MatDiagonalSet_Shell_Private(A,D,1.);
739:     VecDestroy(&d);
740:     if (shell->zrows) {
741:       VecCopy(D,shell->zvals);
742:     }
743:   } else {
744:     MatDiagonalSet_Shell_Private(A,D,1.);
745:     if (shell->zrows) {
746:       VecAXPY(shell->zvals,1.0,D);
747:     }
748:   }
749:   return(0);
750: }

752: PetscErrorCode MatScale_Shell(Mat Y,PetscScalar a)
753: {
754:   Mat_Shell      *shell = (Mat_Shell*)Y->data;

758:   shell->vscale *= a;
759:   shell->vshift *= a;
760:   if (shell->dshift) {
761:     VecScale(shell->dshift,a);
762:   }
763:   shell->axpy_vscale *= a;
764:   if (shell->zrows) {
765:     VecScale(shell->zvals,a);
766:   }
767:   return(0);
768: }

770: static PetscErrorCode MatDiagonalScale_Shell(Mat Y,Vec left,Vec right)
771: {
772:   Mat_Shell      *shell = (Mat_Shell*)Y->data;

776:   if (shell->axpy) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_SUP,"Cannot diagonal scale MATSHELL after MatAXPY operation");
777:   if (left) {
778:     if (!shell->left) {
779:       VecDuplicate(left,&shell->left);
780:       VecCopy(left,shell->left);
781:     } else {
782:       VecPointwiseMult(shell->left,shell->left,left);
783:     }
784:     if (shell->zrows) {
785:       VecPointwiseMult(shell->zvals,shell->zvals,left);
786:     }
787:   }
788:   if (right) {
789:     if (!shell->right) {
790:       VecDuplicate(right,&shell->right);
791:       VecCopy(right,shell->right);
792:     } else {
793:       VecPointwiseMult(shell->right,shell->right,right);
794:     }
795:     if (shell->zrows) {
796:       if (!shell->left_work) {
797:         MatCreateVecs(Y,NULL,&shell->left_work);
798:       }
799:       VecSet(shell->zvals_w,1.0);
800:       VecScatterBegin(shell->zvals_sct_c,right,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD);
801:       VecScatterEnd(shell->zvals_sct_c,right,shell->zvals_w,INSERT_VALUES,SCATTER_FORWARD);
802:       VecPointwiseMult(shell->zvals,shell->zvals,shell->zvals_w);
803:     }
804:   }
805:   return(0);
806: }

808: PetscErrorCode MatAssemblyEnd_Shell(Mat Y,MatAssemblyType t)
809: {
810:   Mat_Shell      *shell = (Mat_Shell*)Y->data;

814:   if (t == MAT_FINAL_ASSEMBLY) {
815:     shell->vshift = 0.0;
816:     shell->vscale = 1.0;
817:     VecDestroy(&shell->dshift);
818:     VecDestroy(&shell->left);
819:     VecDestroy(&shell->right);
820:     MatDestroy(&shell->axpy);
821:     VecScatterDestroy(&shell->zvals_sct_c);
822:     VecScatterDestroy(&shell->zvals_sct_r);
823:     ISDestroy(&shell->zrows);
824:     ISDestroy(&shell->zcols);
825:   }
826:   return(0);
827: }

829: static PetscErrorCode MatMissingDiagonal_Shell(Mat A,PetscBool  *missing,PetscInt *d)
830: {
832:   *missing = PETSC_FALSE;
833:   return(0);
834: }

836: PetscErrorCode MatAXPY_Shell(Mat Y,PetscScalar a,Mat X,MatStructure str)
837: {
838:   Mat_Shell      *shell = (Mat_Shell*)Y->data;

842:   PetscObjectReference((PetscObject)X);
843:   MatDestroy(&shell->axpy);
844:   shell->axpy        = X;
845:   shell->axpy_vscale = a;
846:   return(0);
847: }

849: static struct _MatOps MatOps_Values = {0,
850:                                        0,
851:                                        0,
852:                                        0,
853:                                 /* 4*/ MatMultAdd_Shell,
854:                                        0,
855:                                        MatMultTransposeAdd_Shell,
856:                                        0,
857:                                        0,
858:                                        0,
859:                                 /*10*/ 0,
860:                                        0,
861:                                        0,
862:                                        0,
863:                                        0,
864:                                 /*15*/ 0,
865:                                        0,
866:                                        0,
867:                                        MatDiagonalScale_Shell,
868:                                        0,
869:                                 /*20*/ 0,
870:                                        MatAssemblyEnd_Shell,
871:                                        0,
872:                                        0,
873:                                 /*24*/ MatZeroRows_Shell,
874:                                        0,
875:                                        0,
876:                                        0,
877:                                        0,
878:                                 /*29*/ 0,
879:                                        0,
880:                                        0,
881:                                        0,
882:                                        0,
883:                                 /*34*/ MatDuplicate_Shell,
884:                                        0,
885:                                        0,
886:                                        0,
887:                                        0,
888:                                 /*39*/ MatAXPY_Shell,
889:                                        0,
890:                                        0,
891:                                        0,
892:                                        MatCopy_Shell,
893:                                 /*44*/ 0,
894:                                        MatScale_Shell,
895:                                        MatShift_Shell,
896:                                        MatDiagonalSet_Shell,
897:                                        MatZeroRowsColumns_Shell,
898:                                 /*49*/ 0,
899:                                        0,
900:                                        0,
901:                                        0,
902:                                        0,
903:                                 /*54*/ 0,
904:                                        0,
905:                                        0,
906:                                        0,
907:                                        0,
908:                                 /*59*/ 0,
909:                                        MatDestroy_Shell,
910:                                        0,
911:                                        MatConvertFrom_Shell,
912:                                        0,
913:                                 /*64*/ 0,
914:                                        0,
915:                                        0,
916:                                        0,
917:                                        0,
918:                                 /*69*/ 0,
919:                                        0,
920:                                        MatConvert_Shell,
921:                                        0,
922:                                        0,
923:                                 /*74*/ 0,
924:                                        0,
925:                                        0,
926:                                        0,
927:                                        0,
928:                                 /*79*/ 0,
929:                                        0,
930:                                        0,
931:                                        0,
932:                                        0,
933:                                 /*84*/ 0,
934:                                        0,
935:                                        0,
936:                                        0,
937:                                        0,
938:                                 /*89*/ 0,
939:                                        0,
940:                                        0,
941:                                        0,
942:                                        0,
943:                                 /*94*/ 0,
944:                                        0,
945:                                        0,
946:                                        0,
947:                                        0,
948:                                 /*99*/ 0,
949:                                        0,
950:                                        0,
951:                                        0,
952:                                        0,
953:                                /*104*/ 0,
954:                                        0,
955:                                        0,
956:                                        0,
957:                                        0,
958:                                /*109*/ 0,
959:                                        0,
960:                                        0,
961:                                        0,
962:                                        MatMissingDiagonal_Shell,
963:                                /*114*/ 0,
964:                                        0,
965:                                        0,
966:                                        0,
967:                                        0,
968:                                /*119*/ 0,
969:                                        0,
970:                                        0,
971:                                        0,
972:                                        0,
973:                                /*124*/ 0,
974:                                        0,
975:                                        0,
976:                                        0,
977:                                        0,
978:                                /*129*/ 0,
979:                                        0,
980:                                        0,
981:                                        0,
982:                                        0,
983:                                /*134*/ 0,
984:                                        0,
985:                                        0,
986:                                        0,
987:                                        0,
988:                                /*139*/ 0,
989:                                        0,
990:                                        0
991: };

993: /*MC
994:    MATSHELL - MATSHELL = "shell" - A matrix type to be used to define your own matrix type -- perhaps matrix free.

996:   Level: advanced

998: .seealso: MatCreateShell()
999: M*/

1001: PETSC_EXTERN PetscErrorCode MatCreate_Shell(Mat A)
1002: {
1003:   Mat_Shell      *b;

1007:   PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));

1009:   PetscNewLog(A,&b);
1010:   A->data = (void*)b;

1012:   PetscLayoutSetUp(A->rmap);
1013:   PetscLayoutSetUp(A->cmap);

1015:   b->ctx                 = 0;
1016:   b->vshift              = 0.0;
1017:   b->vscale              = 1.0;
1018:   b->managescalingshifts = PETSC_TRUE;
1019:   A->assembled           = PETSC_TRUE;
1020:   A->preallocated        = PETSC_FALSE;

1022:   PetscObjectChangeTypeName((PetscObject)A,MATSHELL);
1023:   return(0);
1024: }

1026: /*@C
1027:    MatCreateShell - Creates a new matrix class for use with a user-defined
1028:    private data storage format.

1030:   Collective on MPI_Comm

1032:    Input Parameters:
1033: +  comm - MPI communicator
1034: .  m - number of local rows (must be given)
1035: .  n - number of local columns (must be given)
1036: .  M - number of global rows (may be PETSC_DETERMINE)
1037: .  N - number of global columns (may be PETSC_DETERMINE)
1038: -  ctx - pointer to data needed by the shell matrix routines

1040:    Output Parameter:
1041: .  A - the matrix

1043:    Level: advanced

1045:   Usage:
1046: $    extern int mult(Mat,Vec,Vec);
1047: $    MatCreateShell(comm,m,n,M,N,ctx,&mat);
1048: $    MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult);
1049: $    [ Use matrix for operations that have been set ]
1050: $    MatDestroy(mat);

1052:    Notes:
1053:    The shell matrix type is intended to provide a simple class to use
1054:    with KSP (such as, for use with matrix-free methods). You should not
1055:    use the shell type if you plan to define a complete matrix class.

1057:    Fortran Notes:
1058:     To use this from Fortran with a ctx you must write an interface definition for this
1059:     function and for MatShellGetContext() that tells Fortran the Fortran derived data type you are passing
1060:     in as the ctx argument.

1062:    PETSc requires that matrices and vectors being used for certain
1063:    operations are partitioned accordingly.  For example, when
1064:    creating a shell matrix, A, that supports parallel matrix-vector
1065:    products using MatMult(A,x,y) the user should set the number
1066:    of local matrix rows to be the number of local elements of the
1067:    corresponding result vector, y. Note that this is information is
1068:    required for use of the matrix interface routines, even though
1069:    the shell matrix may not actually be physically partitioned.
1070:    For example,

1072: $
1073: $     Vec x, y
1074: $     extern int mult(Mat,Vec,Vec);
1075: $     Mat A
1076: $
1077: $     VecCreateMPI(comm,PETSC_DECIDE,M,&y);
1078: $     VecCreateMPI(comm,PETSC_DECIDE,N,&x);
1079: $     VecGetLocalSize(y,&m);
1080: $     VecGetLocalSize(x,&n);
1081: $     MatCreateShell(comm,m,n,M,N,ctx,&A);
1082: $     MatShellSetOperation(mat,MATOP_MULT,(void(*)(void))mult);
1083: $     MatMult(A,x,y);
1084: $     MatDestroy(&A);
1085: $     VecDestroy(&y);
1086: $     VecDestroy(&x);
1087: $


1090:    MATSHELL handles MatShift(), MatDiagonalSet(), MatDiagonalScale(), MatAXPY(), MatScale(), MatZeroRows() and MatZeroRowsColumns() internally, so these
1091:    operations cannot be overwritten unless MatShellSetManageScalingShifts() is called.


1094:     For rectangular matrices do all the scalings and shifts make sense?

1096:     Developers Notes:
1097:     Regarding shifting and scaling. The general form is

1099:           diag(left)(vscale*A + diag(dshift) + vshift I)diag(right)

1101:       The order you apply the operations is important. For example if you have a dshift then
1102:       apply a MatScale(s) you get s*vscale*A + s*diag(shift). But if you first scale and then shift
1103:       you get s*vscale*A + diag(shift)

1105:           A is the user provided function.

1107:    KSP/PC uses changes in the Mat's "state" to decide if preconditioners need to be rebuilt: PCSetUp() only calls the setup() for
1108:    for the PC implementation if the Mat state has increased from the previous call. Thus to get changes in a MATSHELL to trigger
1109:    an update in the preconditioner you must call MatAssemblyBegin()/MatAssemblyEnd() or PetscObjectStateIncrease((PetscObject)mat);
1110:    each time the MATSHELL matrix has changed.

1112:    Calling MatAssemblyBegin()/MatAssemblyEnd() on a MATSHELL removes any previously supplied shift and scales that were provided
1113:    with MatDiagonalSet(), MatShift(), MatScale(), or MatDiagonalScale().

1115: .keywords: matrix, shell, create

1117: .seealso: MatShellSetOperation(), MatHasOperation(), MatShellGetContext(), MatShellSetContext(), MATSHELL, MatShellSetManageScalingShifts()
1118: @*/
1119: PetscErrorCode  MatCreateShell(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx,Mat *A)
1120: {

1124:   MatCreate(comm,A);
1125:   MatSetSizes(*A,m,n,M,N);
1126:   MatSetType(*A,MATSHELL);
1127:   MatShellSetContext(*A,ctx);
1128:   MatSetUp(*A);
1129:   return(0);
1130: }

1132: /*@
1133:     MatShellSetContext - sets the context for a shell matrix

1135:    Logically Collective on Mat

1137:     Input Parameters:
1138: +   mat - the shell matrix
1139: -   ctx - the context

1141:    Level: advanced

1143:    Fortran Notes:
1144:     To use this from Fortran you must write a Fortran interface definition for this
1145:     function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.

1147: .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation()
1148: @*/
1149: PetscErrorCode  MatShellSetContext(Mat mat,void *ctx)
1150: {
1151:   Mat_Shell      *shell = (Mat_Shell*)mat->data;
1153:   PetscBool      flg;

1157:   PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);
1158:   if (flg) {
1159:     shell->ctx = ctx;
1160:   } else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot attach context to non-shell matrix");
1161:   return(0);
1162: }

1164: /*@
1165:     MatShellSetManageScalingShifts - Allows the user to control the scaling and shift operations of the MATSHELL. Must be called immediately
1166:           after MatCreateShell()

1168:    Logically Collective on Mat

1170:     Input Parameter:
1171: .   mat - the shell matrix

1173:   Level: advanced

1175: .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext(), MatShellSetOperation()
1176: @*/
1177: PetscErrorCode MatShellSetManageScalingShifts(Mat A)
1178: {
1180:   Mat_Shell      *shell;
1181:   PetscBool      flg;

1185:   PetscObjectTypeCompare((PetscObject)A,MATSHELL,&flg);
1186:   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Can only use with MATSHELL matrices");
1187:   shell = (Mat_Shell*)A->data;
1188:   shell->managescalingshifts = PETSC_FALSE;
1189:   A->ops->diagonalset   = NULL;
1190:   A->ops->diagonalscale = NULL;
1191:   A->ops->scale         = NULL;
1192:   A->ops->shift         = NULL;
1193:   A->ops->axpy          = NULL;
1194:   return(0);
1195: }

1197: /*@C
1198:     MatShellTestMult - Compares the multiply routine provided to the MATSHELL with differencing on a given function.

1200:    Logically Collective on Mat

1202:     Input Parameters:
1203: +   mat - the shell matrix
1204: .   f - the function
1205: .   base - differences are computed around this vector, see MatMFFDSetBase(), for Jacobians this is the point at which the Jacobian is being evaluated
1206: -   ctx - an optional context for the function

1208:    Output Parameter:
1209: .   flg - PETSC_TRUE if the multiply is likely correct

1211:    Options Database:
1212: .   -mat_shell_test_mult_view - print if any differences are detected between the products and print the difference

1214:    Level: advanced

1216:    Fortran Notes:
1217:     Not supported from Fortran

1219: .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellTestMultTranspose()
1220: @*/
1221: PetscErrorCode  MatShellTestMult(Mat mat,PetscErrorCode (*f)(void*,Vec,Vec),Vec base,void *ctx,PetscBool *flg)
1222: {
1224:   PetscInt       m,n;
1225:   Mat            mf,Dmf,Dmat,Ddiff;
1226:   PetscReal      Diffnorm,Dmfnorm;
1227:   PetscBool      v = PETSC_FALSE, flag = PETSC_TRUE;

1231:   PetscOptionsHasName(NULL,((PetscObject)mat)->prefix,"-mat_shell_test_mult_view",&v);
1232:   MatGetLocalSize(mat,&m,&n);
1233:   MatCreateMFFD(PetscObjectComm((PetscObject)mat),m,n,PETSC_DECIDE,PETSC_DECIDE,&mf);
1234:   MatMFFDSetFunction(mf,f,ctx);
1235:   MatMFFDSetBase(mf,base,NULL);

1237:   MatComputeOperator(mf,MATAIJ,&Dmf);
1238:   MatComputeOperator(mat,MATAIJ,&Dmat);

1240:   MatDuplicate(Dmat,MAT_COPY_VALUES,&Ddiff);
1241:   MatAXPY(Ddiff,-1.0,Dmf,DIFFERENT_NONZERO_PATTERN);
1242:   MatNorm(Ddiff,NORM_FROBENIUS,&Diffnorm);
1243:   MatNorm(Dmf,NORM_FROBENIUS,&Dmfnorm);
1244:   if (Diffnorm/Dmfnorm > 10*PETSC_SQRT_MACHINE_EPSILON) {
1245:     flag = PETSC_FALSE;
1246:     if (v) {
1247:       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));
1248:       MatViewFromOptions(Ddiff,(PetscObject)mat,"-mat_shell_test_mult_view");
1249:       MatViewFromOptions(Dmf,(PetscObject)mat,"-mat_shell_test_mult_view");
1250:       MatViewFromOptions(Dmat,(PetscObject)mat,"-mat_shell_test_mult_view");
1251:     }
1252:   } else if (v) {
1253:     PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL and matrix free multiple appear to produce the same results\n");
1254:   }
1255:   if (flg) *flg = flag;
1256:   MatDestroy(&Ddiff);
1257:   MatDestroy(&mf);
1258:   MatDestroy(&Dmf);
1259:   MatDestroy(&Dmat);
1260:   return(0);
1261: }

1263: /*@C
1264:     MatShellTestMultTranpose - Compares the multiply transpose routine provided to the MATSHELL with differencing on a given function.

1266:    Logically Collective on Mat

1268:     Input Parameters:
1269: +   mat - the shell matrix
1270: .   f - the function
1271: .   base - differences are computed around this vector, see MatMFFDSetBase(), for Jacobians this is the point at which the Jacobian is being evaluated
1272: -   ctx - an optional context for the function

1274:    Output Parameter:
1275: .   flg - PETSC_TRUE if the multiply is likely correct

1277:    Options Database:
1278: .   -mat_shell_test_mult_view - print if any differences are detected between the products and print the difference

1280:    Level: advanced

1282:    Fortran Notes:
1283:     Not supported from Fortran

1285: .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellTestMult()
1286: @*/
1287: PetscErrorCode  MatShellTestMultTranspose(Mat mat,PetscErrorCode (*f)(void*,Vec,Vec),Vec base,void *ctx,PetscBool *flg)
1288: {
1290:   Vec            x,y,z;
1291:   PetscInt       m,n,M,N;
1292:   Mat            mf,Dmf,Dmat,Ddiff;
1293:   PetscReal      Diffnorm,Dmfnorm;
1294:   PetscBool      v = PETSC_FALSE, flag = PETSC_TRUE;

1298:   PetscOptionsHasName(NULL,((PetscObject)mat)->prefix,"-mat_shell_test_mult_transpose_view",&v);
1299:   MatCreateVecs(mat,&x,&y);
1300:   VecDuplicate(y,&z);
1301:   MatGetLocalSize(mat,&m,&n);
1302:   MatGetSize(mat,&M,&N);
1303:   MatCreateMFFD(PetscObjectComm((PetscObject)mat),m,n,M,N,&mf);
1304:   MatMFFDSetFunction(mf,f,ctx);
1305:   MatMFFDSetBase(mf,base,NULL);
1306:   MatComputeOperator(mf,MATAIJ,&Dmf);
1307:   MatTranspose(Dmf,MAT_INPLACE_MATRIX,&Dmf);
1308:   MatComputeOperatorTranspose(mat,MATAIJ,&Dmat);

1310:   MatDuplicate(Dmat,MAT_COPY_VALUES,&Ddiff);
1311:   MatAXPY(Ddiff,-1.0,Dmf,DIFFERENT_NONZERO_PATTERN);
1312:   MatNorm(Ddiff,NORM_FROBENIUS,&Diffnorm);
1313:   MatNorm(Dmf,NORM_FROBENIUS,&Dmfnorm);
1314:   if (Diffnorm/Dmfnorm > 10*PETSC_SQRT_MACHINE_EPSILON) {
1315:     flag = PETSC_FALSE;
1316:     if (v) {
1317:       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));
1318:       MatViewFromOptions(Ddiff,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");
1319:       MatViewFromOptions(Dmf,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");
1320:       MatViewFromOptions(Dmat,(PetscObject)mat,"-mat_shell_test_mult_transpose_view");
1321:     }
1322:   } else if (v) {
1323:     PetscPrintf(PetscObjectComm((PetscObject)mat),"MATSHELL transpose and matrix free multiple appear to produce the same results\n");
1324:   }
1325:   if (flg) *flg = flag;
1326:   MatDestroy(&mf);
1327:   MatDestroy(&Dmat);
1328:   MatDestroy(&Ddiff);
1329:   MatDestroy(&Dmf);
1330:   VecDestroy(&x);
1331:   VecDestroy(&y);
1332:   VecDestroy(&z);
1333:   return(0);
1334: }

1336: /*@C
1337:     MatShellSetOperation - Allows user to set a matrix operation for a shell matrix.

1339:    Logically Collective on Mat

1341:     Input Parameters:
1342: +   mat - the shell matrix
1343: .   op - the name of the operation
1344: -   f - the function that provides the operation.

1346:    Level: advanced

1348:     Usage:
1349: $      extern PetscErrorCode usermult(Mat,Vec,Vec);
1350: $      MatCreateShell(comm,m,n,M,N,ctx,&A);
1351: $      MatShellSetOperation(A,MATOP_MULT,(void(*)(void))usermult);

1353:     Notes:
1354:     See the file include/petscmat.h for a complete list of matrix
1355:     operations, which all have the form MATOP_<OPERATION>, where
1356:     <OPERATION> is the name (in all capital letters) of the
1357:     user interface routine (e.g., MatMult() -> MATOP_MULT).

1359:     All user-provided functions (except for MATOP_DESTROY) should have the same calling
1360:     sequence as the usual matrix interface routines, since they
1361:     are intended to be accessed via the usual matrix interface
1362:     routines, e.g.,
1363: $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)

1365:     In particular each function MUST return an error code of 0 on success and
1366:     nonzero on failure.

1368:     Within each user-defined routine, the user should call
1369:     MatShellGetContext() to obtain the user-defined context that was
1370:     set by MatCreateShell().

1372:     Fortran Notes:
1373:     For MatCreateVecs() the user code should check if the input left or right matrix is -1 and in that case not
1374:        generate a matrix. See src/mat/examples/tests/ex120f.F

1376:     Use MatSetOperation() to set an operation for any matrix type

1378: .keywords: matrix, shell, set, operation

1380: .seealso: MatCreateShell(), MatShellGetContext(), MatShellGetOperation(), MatShellSetContext(), MatSetOperation(), MatShellSetManageScalingShifts()
1381: @*/
1382: PetscErrorCode MatShellSetOperation(Mat mat,MatOperation op,void (*f)(void))
1383: {
1384:   PetscBool      flg;
1385:   Mat_Shell      *shell;

1390:   PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);
1391:   if (!flg) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Can only use with MATSHELL matrices");
1392:   shell = (Mat_Shell*)mat->data;

1394:   switch (op) {
1395:   case MATOP_DESTROY:
1396:     shell->ops->destroy = (PetscErrorCode (*)(Mat))f;
1397:     break;
1398:   case MATOP_VIEW:
1399:     if (!mat->ops->viewnative) {
1400:       mat->ops->viewnative = mat->ops->view;
1401:     }
1402:     mat->ops->view = (PetscErrorCode (*)(Mat,PetscViewer))f;
1403:     break;
1404:   case MATOP_COPY:
1405:     shell->ops->copy = (PetscErrorCode (*)(Mat,Mat,MatStructure))f;
1406:     break;
1407:   case MATOP_DIAGONAL_SET:
1408:   case MATOP_DIAGONAL_SCALE:
1409:   case MATOP_SHIFT:
1410:   case MATOP_SCALE:
1411:   case MATOP_AXPY:
1412:   case MATOP_ZERO_ROWS:
1413:   case MATOP_ZERO_ROWS_COLUMNS:
1414:     if (shell->managescalingshifts) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"MATSHELL is managing scalings and shifts, see MatShellSetManageScalingShifts()");
1415:     (((void(**)(void))mat->ops)[op]) = f;
1416:     break;
1417:   case MATOP_GET_DIAGONAL:
1418:     if (shell->managescalingshifts) {
1419:       shell->ops->getdiagonal = (PetscErrorCode (*)(Mat,Vec))f;
1420:       mat->ops->getdiagonal   = MatGetDiagonal_Shell;
1421:     } else {
1422:       shell->ops->getdiagonal = NULL;
1423:       mat->ops->getdiagonal   = (PetscErrorCode (*)(Mat,Vec))f;
1424:     }
1425:     break;
1426:   case MATOP_MULT:
1427:     if (shell->managescalingshifts) {
1428:       shell->ops->mult = (PetscErrorCode (*)(Mat,Vec,Vec))f;
1429:       mat->ops->mult   = MatMult_Shell;
1430:     } else {
1431:       shell->ops->mult = NULL;
1432:       mat->ops->mult   = (PetscErrorCode (*)(Mat,Vec,Vec))f;
1433:     }
1434:     break;
1435:   case MATOP_MULT_TRANSPOSE:
1436:     if (shell->managescalingshifts) {
1437:       shell->ops->multtranspose = (PetscErrorCode (*)(Mat,Vec,Vec))f;
1438:       mat->ops->multtranspose   = MatMultTranspose_Shell;
1439:     } else {
1440:       shell->ops->multtranspose = NULL;
1441:       mat->ops->multtranspose   = (PetscErrorCode (*)(Mat,Vec,Vec))f;
1442:     }
1443:     break;
1444:   default:
1445:     (((void(**)(void))mat->ops)[op]) = f;
1446:     break;
1447:   }
1448:   return(0);
1449: }

1451: /*@C
1452:     MatShellGetOperation - Gets a matrix function for a shell matrix.

1454:     Not Collective

1456:     Input Parameters:
1457: +   mat - the shell matrix
1458: -   op - the name of the operation

1460:     Output Parameter:
1461: .   f - the function that provides the operation.

1463:     Level: advanced

1465:     Notes:
1466:     See the file include/petscmat.h for a complete list of matrix
1467:     operations, which all have the form MATOP_<OPERATION>, where
1468:     <OPERATION> is the name (in all capital letters) of the
1469:     user interface routine (e.g., MatMult() -> MATOP_MULT).

1471:     All user-provided functions have the same calling
1472:     sequence as the usual matrix interface routines, since they
1473:     are intended to be accessed via the usual matrix interface
1474:     routines, e.g.,
1475: $       MatMult(Mat,Vec,Vec) -> usermult(Mat,Vec,Vec)

1477:     Within each user-defined routine, the user should call
1478:     MatShellGetContext() to obtain the user-defined context that was
1479:     set by MatCreateShell().

1481: .keywords: matrix, shell, set, operation

1483: .seealso: MatCreateShell(), MatShellGetContext(), MatShellSetOperation(), MatShellSetContext()
1484: @*/
1485: PetscErrorCode MatShellGetOperation(Mat mat,MatOperation op,void(**f)(void))
1486: {
1487:   PetscBool      flg;
1488:   Mat_Shell      *shell;

1493:   PetscObjectTypeCompare((PetscObject)mat,MATSHELL,&flg);
1494:   if (!flg) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Can only use with MATSHELL matrices");
1495:   shell = (Mat_Shell*)mat->data;

1497:   switch (op) {
1498:   case MATOP_DESTROY:
1499:     *f = (void (*)(void))shell->ops->destroy;
1500:     break;
1501:   case MATOP_VIEW:
1502:     *f = (void (*)(void))mat->ops->view;
1503:     break;
1504:   case MATOP_COPY:
1505:     *f = (void (*)(void))shell->ops->copy;
1506:     break;
1507:   case MATOP_DIAGONAL_SET:
1508:   case MATOP_DIAGONAL_SCALE:
1509:   case MATOP_SHIFT:
1510:   case MATOP_SCALE:
1511:   case MATOP_AXPY:
1512:   case MATOP_ZERO_ROWS:
1513:   case MATOP_ZERO_ROWS_COLUMNS:
1514:     *f = (((void (**)(void))mat->ops)[op]);
1515:     break;
1516:   case MATOP_GET_DIAGONAL:
1517:     if (shell->ops->getdiagonal)
1518:       *f = (void (*)(void))shell->ops->getdiagonal;
1519:     else
1520:       *f = (((void (**)(void))mat->ops)[op]);
1521:     break;
1522:   case MATOP_MULT:
1523:     if (shell->ops->mult)
1524:       *f = (void (*)(void))shell->ops->mult;
1525:     else
1526:       *f = (((void (**)(void))mat->ops)[op]);
1527:     break;
1528:   case MATOP_MULT_TRANSPOSE:
1529:     if (shell->ops->multtranspose)
1530:       *f = (void (*)(void))shell->ops->multtranspose;
1531:     else
1532:       *f = (((void (**)(void))mat->ops)[op]);
1533:     break;
1534:   default:
1535:     *f = (((void (**)(void))mat->ops)[op]);
1536:   }
1537:   return(0);
1538: }