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.

  6:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  7:       SLEPc - Scalable Library for Eigenvalue Problem Computations
  8:       Copyright (c) 2002-2007, Universidad Politecnica de Valencia, Spain

 10:       This file is part of SLEPc. See the README file for conditions of use
 11:       and additional information.
 12:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 13: */


 16:  #include src/st/stimpl.h

 20: PetscErrorCode STMatShellMult(Mat A,Vec x,Vec y)
 21: {
 23:   ST             ctx;

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

 28:   MatMult(ctx->A,x,y);
 29:   if (ctx->sigma != 0.0) {
 30:     if (ctx->B) {  /* y = (A - sB) x */
 31:       MatMult(ctx->B,x,ctx->w);
 32:       VecAXPY(y,-ctx->sigma,ctx->w);
 33:     } else {    /* y = (A - sI) x */
 34:       VecAXPY(y,-ctx->sigma,x);
 35:     }
 36:   }
 37:   return(0);
 38: }

 42: PetscErrorCode STMatShellMultTranspose(Mat A,Vec x,Vec y)
 43: {
 45:   ST             ctx;

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

 50:   MatMultTranspose(ctx->A,x,y);
 51:   if (ctx->sigma != 0.0) {
 52:     if (ctx->B) {  /* y = (A - sB) x */
 53:       MatMultTranspose(ctx->B,x,ctx->w);
 54:       VecAXPY(y,-ctx->sigma,ctx->w);
 55:     } else {    /* y = (A - sI) x */
 56:       VecAXPY(y,-ctx->sigma,x);
 57:     }
 58:   }
 59:   return(0);
 60: }

 64: PetscErrorCode STMatShellGetDiagonal(Mat A,Vec diag)
 65: {
 67:   ST             ctx;
 68:   Vec            diagb;

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

 73:   MatGetDiagonal(ctx->A,diag);
 74:   if (ctx->sigma != 0.0) {
 75:     if (ctx->B) {
 76:       VecDuplicate(diag,&diagb);
 77:       MatGetDiagonal(ctx->B,diagb);
 78:       VecAXPY(diag,-ctx->sigma,diagb);
 79:       VecDestroy(diagb);
 80:     } else {
 81:       VecShift(diag,-ctx->sigma);
 82:     }
 83:   }
 84:   return(0);
 85: }

 89: PetscErrorCode STMatShellCreate(ST st,Mat *mat)
 90: {
 92:   PetscInt       n, m, N, M;
 93:   PetscTruth     hasA, hasB;

 96:   MatGetSize(st->A,&M,&N);
 97:   MatGetLocalSize(st->A,&m,&n);
 98:   MatCreateShell(st->comm,m,n,M,N,(void*)st,mat);
 99:   MatShellSetOperation(*mat,MATOP_MULT,(void(*)(void))STMatShellMult);
100:   MatShellSetOperation(*mat,MATOP_MULT_TRANSPOSE,(void(*)(void))STMatShellMultTranspose);

102:   MatHasOperation(st->A,MATOP_GET_DIAGONAL,&hasA);
103:   if (st->B) { MatHasOperation(st->B,MATOP_GET_DIAGONAL,&hasB); }
104:   if ( (hasA && !st->B) || (hasA && hasB) ) {
105:     MatShellSetOperation(*mat,MATOP_GET_DIAGONAL,(void(*)(void))STMatShellGetDiagonal);
106:   }

108:   return(0);
109: }