Actual source code: test1.c

  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2012, Universitat Politecnica de Valencia, Spain

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

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

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

 22: static char help[] = "Test ST with shell matrices as problem matrices.\n\n";

 24: #include <slepceps.h>

 26: static PetscErrorCode MatShift_Shell(Mat S,PetscScalar shift);
 27: static PetscErrorCode MatGetDiagonal_Shell(Mat S,Vec diag);
 28: static PetscErrorCode MatMultTranspose_Shell(Mat S,Vec x,Vec y);
 29: static PetscErrorCode MatMult_Shell(Mat S,Vec x,Vec y);

 33: int main(int argc,char **argv)
 34: {
 35:   Mat            A, S;        /* problem matrix */
 36:   EPS            eps;         /* eigenproblem solver context */
 37:   const EPSType  type;
 38:   PetscScalar    value[3];
 39:   PetscInt       n=30,i,Istart,Iend,col[3],nev;
 40:   PetscBool      FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;

 43:   SlepcInitialize(&argc,&argv,(char*)0,help);

 45:   PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
 46:   PetscPrintf(PETSC_COMM_WORLD,"\n1-D Laplacian Eigenproblem, n=%D\n\n",n);

 48:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 49:      Compute the operator matrix that defines the eigensystem, Ax=kx
 50:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 52:   MatCreate(PETSC_COMM_WORLD,&A);
 53:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
 54:   MatSetFromOptions(A);
 55:   MatSetUp(A);
 56: 
 57:   MatGetOwnershipRange(A,&Istart,&Iend);
 58:   if (Istart==0) FirstBlock=PETSC_TRUE;
 59:   if (Iend==n) LastBlock=PETSC_TRUE;
 60:   value[0]=-1.0; value[1]=2.0; value[2]=-1.0;
 61:   for (i=(FirstBlock? Istart+1: Istart); i<(LastBlock? Iend-1: Iend); i++) {
 62:     col[0]=i-1; col[1]=i; col[2]=i+1;
 63:     MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
 64:   }
 65:   if (LastBlock) {
 66:     i=n-1; col[0]=n-2; col[1]=n-1;
 67:     MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
 68:   }
 69:   if (FirstBlock) {
 70:     i=0; col[0]=0; col[1]=1; value[0]=2.0; value[1]=-1.0;
 71:     MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
 72:   }

 74:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 75:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 77:   /* Create the shell version of A */
 78:   MatCreateShell(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,&A,&S);
 79:   MatSetFromOptions(S);
 80:   MatShellSetOperation(S,MATOP_MULT,(void(*)())MatMult_Shell);
 81:   MatShellSetOperation(S,MATOP_MULT_TRANSPOSE,(void(*)())MatMultTranspose_Shell);
 82:   MatShellSetOperation(S,MATOP_GET_DIAGONAL,(void(*)())MatGetDiagonal_Shell);
 83:   MatShellSetOperation(S,MATOP_SHIFT,(void(*)())MatShift_Shell);

 85:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 86:                 Create the eigensolver and set various options
 87:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 89:   EPSCreate(PETSC_COMM_WORLD,&eps);
 90:   EPSSetOperators(eps,S,PETSC_NULL);
 91:   EPSSetProblemType(eps,EPS_HEP);
 92:   EPSSetFromOptions(eps);

 94:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 95:                       Solve the eigensystem
 96:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 98:   EPSSolve(eps);
 99:   EPSGetType(eps,&type);
100:   PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
101:   EPSGetDimensions(eps,&nev,PETSC_NULL,PETSC_NULL);
102:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);

104:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
105:                     Display solution and clean up
106:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

108:   EPSPrintSolution(eps,PETSC_NULL);
109:   EPSDestroy(&eps);
110:   MatDestroy(&A);
111:   MatDestroy(&S);
112:   SlepcFinalize();
113:   return 0;
114: }

118: static PetscErrorCode MatMult_Shell(Mat S,Vec x,Vec y)
119: {
120:   PetscErrorCode    ierr;
121:   Mat               *A;
122: 
124:   MatShellGetContext(S,&A);
125:   MatMult(*A,x,y);
126:   return(0);
127: }

131: static PetscErrorCode MatMultTranspose_Shell(Mat S,Vec x,Vec y)
132: {
133:   PetscErrorCode    ierr;
134:   Mat               *A;
135: 
137:   MatShellGetContext(S,&A);
138:   MatMultTranspose(*A,x,y);
139:   return(0);
140: }

144: static PetscErrorCode MatGetDiagonal_Shell(Mat S,Vec diag)
145: {
146:   PetscErrorCode    ierr;
147:   Mat               *A;
148: 
150:   MatShellGetContext(S,&A);
151:   MatGetDiagonal(*A,diag);
152:   return(0);
153: }

157: static PetscErrorCode MatShift_Shell(Mat S,PetscScalar shift)
158: {
159:   PetscErrorCode    ierr;
160:   Mat               *A;
161: 
163:   MatShellGetContext(S,&A);
164:   MatShift(*A,shift);
165:   return(0);
166: }
167: