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-2011, Universitat 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: static PetscErrorCode STMatShellMult(Mat A,Vec x,Vec y)
 32: {
 34:   ST             ctx;

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

 52: static PetscErrorCode STMatShellMultTranspose(Mat A,Vec x,Vec y)
 53: {
 55:   ST             ctx;

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

 73: static PetscErrorCode STMatShellGetDiagonal(Mat A,Vec diag)
 74: {
 76:   ST             ctx;
 77:   Vec            diagb;

 80:   MatShellGetContext(A,(void**)&ctx);
 81:   MatGetDiagonal(ctx->A,diag);
 82:   if (ctx->sigma != 0.0) {
 83:     if (ctx->B) {
 84:       VecDuplicate(diag,&diagb);
 85:       MatGetDiagonal(ctx->B,diagb);
 86:       VecAXPY(diag,-ctx->sigma,diagb);
 87:       VecDestroy(&diagb);
 88:     } else {
 89:       VecShift(diag,-ctx->sigma);
 90:     }
 91:   }
 92:   return(0);
 93: }

 97: PetscErrorCode STMatShellCreate(ST st,Mat *mat)
 98: {
100:   PetscInt       n,m,N,M;
101:   PetscBool      hasA,hasB;

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

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