Actual source code: test3.c

slepc-3.5.2 2014-10-10
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2014, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.

  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 two matrices.\n\n";

 24: #include <slepcst.h>

 28: int main(int argc,char **argv)
 29: {
 30:   Mat            A,B,M,mat[2];
 31:   ST             st;
 32:   Vec            v,w;
 33:   STType         type;
 34:   PetscScalar    value[3],sigma,tau;
 35:   PetscInt       n=10,i,Istart,Iend,col[3];
 36:   PetscBool      FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;

 39:   SlepcInitialize(&argc,&argv,(char*)0,help);
 40:   PetscOptionsGetInt(NULL,"-n",&n,NULL);
 41:   PetscPrintf(PETSC_COMM_WORLD,"\n1-D Laplacian plus diagonal, n=%D\n\n",n);

 43:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 44:      Compute the operator matrix for the 1-D Laplacian
 45:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 47:   MatCreate(PETSC_COMM_WORLD,&A);
 48:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
 49:   MatSetFromOptions(A);
 50:   MatSetUp(A);

 52:   MatCreate(PETSC_COMM_WORLD,&B);
 53:   MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n);
 54:   MatSetFromOptions(B);
 55:   MatSetUp(B);

 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:     MatSetValue(B,i,i,(PetscScalar)i,INSERT_VALUES);
 65:   }
 66:   if (LastBlock) {
 67:     i=n-1; col[0]=n-2; col[1]=n-1;
 68:     MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
 69:     MatSetValue(B,i,i,(PetscScalar)i,INSERT_VALUES);
 70:   }
 71:   if (FirstBlock) {
 72:     i=0; col[0]=0; col[1]=1; value[0]=2.0; value[1]=-1.0;
 73:     MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
 74:     MatSetValue(B,i,i,-1.0,INSERT_VALUES);
 75:   }

 77:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 78:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 79:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 80:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
 81:   MatGetVecs(A,&v,&w);
 82:   VecSet(v,1.0);

 84:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 85:                 Create the spectral transformation object
 86:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 88:   STCreate(PETSC_COMM_WORLD,&st);
 89:   mat[0] = A;
 90:   mat[1] = B;
 91:   STSetOperators(st,2,mat);
 92:   STSetFromOptions(st);

 94:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 95:               Apply the transformed operator for several ST's
 96:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 98:   /* shift, sigma=0.0 */
 99:   STSetUp(st);
100:   STGetType(st,&type);
101:   PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type);
102:   STApply(st,v,w);
103:   VecView(w,NULL);

105:   /* shift, sigma=0.1 */
106:   sigma = 0.1;
107:   STSetShift(st,sigma);
108:   STGetShift(st,&sigma);
109:   PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma));
110:   STApply(st,v,w);
111:   VecView(w,NULL);

113:   /* sinvert, sigma=0.1 */
114:   STPostSolve(st);   /* undo changes if inplace */
115:   STSetType(st,STSINVERT);
116:   STGetType(st,&type);
117:   PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type);
118:   STGetShift(st,&sigma);
119:   PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma));
120:   STApply(st,v,w);
121:   VecView(w,NULL);

123:   /* sinvert, sigma=-0.5 */
124:   sigma = -0.5;
125:   STSetShift(st,sigma);
126:   STGetShift(st,&sigma);
127:   PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma));
128:   STApply(st,v,w);
129:   VecView(w,NULL);

131:   /* cayley, sigma=-0.5, tau=-0.5 (equal to sigma by default) */
132:   STPostSolve(st);   /* undo changes if inplace */
133:   STSetType(st,STCAYLEY);
134:   STSetUp(st);
135:   STGetType(st,&type);
136:   PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type);
137:   STGetShift(st,&sigma);
138:   STCayleyGetAntishift(st,&tau);
139:   PetscPrintf(PETSC_COMM_WORLD,"With shift=%g, antishift=%g\n",(double)PetscRealPart(sigma),(double)PetscRealPart(tau));
140:   STApply(st,v,w);
141:   VecView(w,NULL);

143:   /* cayley, sigma=1.1, tau=1.1 (still equal to sigma) */
144:   sigma = 1.1;
145:   STSetShift(st,sigma);
146:   STGetShift(st,&sigma);
147:   STCayleyGetAntishift(st,&tau);
148:   PetscPrintf(PETSC_COMM_WORLD,"With shift=%g, antishift=%g\n",(double)PetscRealPart(sigma),(double)PetscRealPart(tau));
149:   STApply(st,v,w);
150:   VecView(w,NULL);

152:   /* cayley, sigma=1.1, tau=-1.0 */
153:   tau = -1.0;
154:   STCayleySetAntishift(st,tau);
155:   STGetShift(st,&sigma);
156:   STCayleyGetAntishift(st,&tau);
157:   PetscPrintf(PETSC_COMM_WORLD,"With shift=%g, antishift=%g\n",(double)PetscRealPart(sigma),(double)PetscRealPart(tau));
158:   STApply(st,v,w);
159:   VecView(w,NULL);

161:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
162:                   Check inner product matrix in Cayley
163:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
164:   STGetBilinearForm(st,&M);
165:   MatMult(M,v,w);
166:   VecView(w,NULL);

168:   STDestroy(&st);
169:   MatDestroy(&A);
170:   MatDestroy(&B);
171:   MatDestroy(&M);
172:   VecDestroy(&v);
173:   VecDestroy(&w);
174:   SlepcFinalize();
175:   return 0;
176: }