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-2009, Universidad Politecnica de Valencia, Spain

 10:    This file is part of SLEPc.
 11:       
 12:    SLEPc is free software: you can redistribute it and/or modify it under  the
 13:    terms of version 3 of the GNU Lesser General Public License as published by
 14:    the Free Software Foundation.

 16:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY 
 17:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS 
 18:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for 
 19:    more details.

 21:    You  should have received a copy of the GNU Lesser General  Public  License
 22:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 23:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 24: */


 27:  #include private/stimpl.h

 31: PetscErrorCode STMatShellMult(Mat A,Vec x,Vec y)
 32: {
 34:   ST             ctx;

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

 39:   MatMult(ctx->A,x,y);
 40:   if (ctx->sigma != 0.0) {
 41:     if (ctx->B) {  /* y = (A - sB) x */
 42:       MatMult(ctx->B,x,ctx->w);
 43:       VecAXPY(y,-ctx->sigma,ctx->w);
 44:     } else {    /* y = (A - sI) x */
 45:       VecAXPY(y,-ctx->sigma,x);
 46:     }
 47:   }
 48:   return(0);
 49: }

 53: PetscErrorCode STMatShellMultTranspose(Mat A,Vec x,Vec y)
 54: {
 56:   ST             ctx;

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

 61:   MatMultTranspose(ctx->A,x,y);
 62:   if (ctx->sigma != 0.0) {
 63:     if (ctx->B) {  /* y = (A - sB) x */
 64:       MatMultTranspose(ctx->B,x,ctx->w);
 65:       VecAXPY(y,-ctx->sigma,ctx->w);
 66:     } else {    /* y = (A - sI) x */
 67:       VecAXPY(y,-ctx->sigma,x);
 68:     }
 69:   }
 70:   return(0);
 71: }

 75: PetscErrorCode STMatShellGetDiagonal(Mat A,Vec diag)
 76: {
 78:   ST             ctx;
 79:   Vec            diagb;

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

 84:   MatGetDiagonal(ctx->A,diag);
 85:   if (ctx->sigma != 0.0) {
 86:     if (ctx->B) {
 87:       VecDuplicate(diag,&diagb);
 88:       MatGetDiagonal(ctx->B,diagb);
 89:       VecAXPY(diag,-ctx->sigma,diagb);
 90:       VecDestroy(diagb);
 91:     } else {
 92:       VecShift(diag,-ctx->sigma);
 93:     }
 94:   }
 95:   return(0);
 96: }

100: PetscErrorCode STMatShellCreate(ST st,Mat *mat)
101: {
103:   PetscInt       n, m, N, M;
104:   PetscTruth     hasA, hasB;

107:   MatGetSize(st->A,&M,&N);
108:   MatGetLocalSize(st->A,&m,&n);
109:   MatCreateShell(((PetscObject)st)->comm,m,n,M,N,(void*)st,mat);
110:   MatShellSetOperation(*mat,MATOP_MULT,(void(*)(void))STMatShellMult);
111:   MatShellSetOperation(*mat,MATOP_MULT_TRANSPOSE,(void(*)(void))STMatShellMultTranspose);

113:   MatHasOperation(st->A,MATOP_GET_DIAGONAL,&hasA);
114:   if (st->B) { MatHasOperation(st->B,MATOP_GET_DIAGONAL,&hasB); }
115:   if ( (hasA && !st->B) || (hasA && hasB) ) {
116:     MatShellSetOperation(*mat,MATOP_GET_DIAGONAL,(void(*)(void))STMatShellGetDiagonal);
117:   }

119:   return(0);
120: }