Actual source code: test3.c

slepc-3.7.3 2016-09-29
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2016, 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    sigma,tau;
 35:   PetscInt       n=10,i,Istart,Iend;

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

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

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

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

 56:   MatGetOwnershipRange(A,&Istart,&Iend);
 57:   for (i=Istart;i<Iend;i++) {
 58:     MatSetValue(A,i,i,2.0,INSERT_VALUES);
 59:     if (i>0) {
 60:       MatSetValue(A,i,i-1,-1.0,INSERT_VALUES);
 61:       MatSetValue(B,i,i,(PetscScalar)i,INSERT_VALUES);
 62:     } else {
 63:       MatSetValue(B,i,i,-1.0,INSERT_VALUES);
 64:     }
 65:     if (i<n-1) { MatSetValue(A,i,i+1,-1.0,INSERT_VALUES); }
 66:   }
 67:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 68:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 69:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 70:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
 71:   MatCreateVecs(A,&v,&w);
 72:   VecSet(v,1.0);

 74:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 75:                 Create the spectral transformation object
 76:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 78:   STCreate(PETSC_COMM_WORLD,&st);
 79:   mat[0] = A;
 80:   mat[1] = B;
 81:   STSetOperators(st,2,mat);
 82:   STSetTransform(st,PETSC_TRUE);
 83:   STSetFromOptions(st);

 85:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 86:               Apply the transformed operator for several ST's
 87:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 89:   /* shift, sigma=0.0 */
 90:   STSetUp(st);
 91:   STGetType(st,&type);
 92:   PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type);
 93:   STApply(st,v,w);
 94:   VecView(w,NULL);

 96:   /* shift, sigma=0.1 */
 97:   sigma = 0.1;
 98:   STSetShift(st,sigma);
 99:   STGetShift(st,&sigma);
100:   PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma));
101:   STApply(st,v,w);
102:   VecView(w,NULL);

104:   /* sinvert, sigma=0.1 */
105:   STPostSolve(st);   /* undo changes if inplace */
106:   STSetType(st,STSINVERT);
107:   STGetType(st,&type);
108:   PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type);
109:   STGetShift(st,&sigma);
110:   PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma));
111:   STApply(st,v,w);
112:   VecView(w,NULL);

114:   /* sinvert, sigma=-0.5 */
115:   sigma = -0.5;
116:   STSetShift(st,sigma);
117:   STGetShift(st,&sigma);
118:   PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma));
119:   STApply(st,v,w);
120:   VecView(w,NULL);

122:   /* cayley, sigma=-0.5, tau=-0.5 (equal to sigma by default) */
123:   STPostSolve(st);   /* undo changes if inplace */
124:   STSetType(st,STCAYLEY);
125:   STSetUp(st);
126:   STGetType(st,&type);
127:   PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type);
128:   STGetShift(st,&sigma);
129:   STCayleyGetAntishift(st,&tau);
130:   PetscPrintf(PETSC_COMM_WORLD,"With shift=%g, antishift=%g\n",(double)PetscRealPart(sigma),(double)PetscRealPart(tau));
131:   STApply(st,v,w);
132:   VecView(w,NULL);

134:   /* cayley, sigma=1.1, tau=1.1 (still equal to sigma) */
135:   sigma = 1.1;
136:   STSetShift(st,sigma);
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.0 */
144:   tau = -1.0;
145:   STCayleySetAntishift(st,tau);
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:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
153:                   Check inner product matrix in Cayley
154:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
155:   STGetBilinearForm(st,&M);
156:   MatMult(M,v,w);
157:   VecView(w,NULL);

159:   STDestroy(&st);
160:   MatDestroy(&A);
161:   MatDestroy(&B);
162:   MatDestroy(&M);
163:   VecDestroy(&v);
164:   VecDestroy(&w);
165:   SlepcFinalize();
166:   return ierr;
167: }