Actual source code: shellmat.c

  1: /*
  2:       This file contains the subroutines which implement various operations
  3:       of the matrix associated to the shift-and-invert technique for eigenvalue
  4:       problems, and also a subroutine to create it.
  5: */

 7:  #include src/st/stimpl.h

 11: PetscErrorCode STMatShellMult(Mat A,Vec x,Vec y)
 12: {
 13:   PetscErrorCode ierr;
 14:   ST             ctx;
 15:   PetscScalar    alpha;

 18:   MatShellGetContext(A,(void**)&ctx);

 20:   MatMult(ctx->A,x,y);
 21:   if (ctx->sigma != 0.0) {
 22:     alpha = -ctx->sigma;
 23:     if (ctx->B) {  /* y = (A - sB) x */
 24:       MatMult(ctx->B,x,ctx->w);
 25:       VecAXPY(&alpha,ctx->w,y);
 26:     }
 27:     else {    /* y = (A - sI) x */
 28:     VecAXPY(&alpha,x,y);
 29:     }
 30:   }
 31:   return(0);
 32: }

 36: PetscErrorCode STMatShellGetDiagonal(Mat A,Vec diag)
 37: {
 38:   PetscErrorCode ierr;
 39:   ST             ctx;
 40:   PetscScalar    alpha;
 41:   Vec            diagb;

 44:   MatShellGetContext(A,(void**)&ctx);

 46:   MatGetDiagonal(ctx->A,diag);
 47:   if (ctx->sigma != 0.0) {
 48:     alpha = -ctx->sigma;
 49:     if (ctx->B) {
 50:       VecDuplicate(diag,&diagb);
 51:       MatGetDiagonal(ctx->B,diagb);
 52:       VecAXPY(&alpha,diagb,diag);
 53:       VecDestroy(diagb);
 54:     }
 55:     else {
 56:       VecShift(&alpha,diag);
 57:     }
 58:   }
 59:   return(0);
 60: }

 64: PetscErrorCode STMatShellCreate(ST st,Mat *mat)
 65: {
 66:   PetscErrorCode ierr;
 67:   int            n, m, N, M;
 68:   PetscTruth     hasA, hasB;

 71:   MatGetSize(st->A,&M,&N);
 72:   MatGetLocalSize(st->A,&m,&n);
 73:   MatCreateShell(st->comm,m,n,M,N,(void*)st,mat);
 74:   MatShellSetOperation(*mat,MATOP_MULT,(void(*)(void))STMatShellMult);

 76:   MatHasOperation(st->A,MATOP_GET_DIAGONAL,&hasA);
 77:   if (st->B) { MatHasOperation(st->B,MATOP_GET_DIAGONAL,&hasB); }
 78:   if ( (hasA && !st->B) || (hasA && hasB) ) {
 79:     MatShellSetOperation(*mat,MATOP_GET_DIAGONAL,(void(*)(void))STMatShellGetDiagonal);
 80:   }

 82:   return(0);
 83: }