Actual source code: test2.c

slepc-3.7.2 2016-07-19
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 one matrix.\n\n";

 24: #include <slepcst.h>

 28: int main(int argc,char **argv)
 29: {
 30:   Mat            A,B,mat[1];
 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, 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:   MatGetOwnershipRange(A,&Istart,&Iend);
 52:   for (i=Istart;i<Iend;i++) {
 53:     if (i>0) { MatSetValue(A,i,i-1,-1.0,INSERT_VALUES); }
 54:     if (i<n-1) { MatSetValue(A,i,i+1,-1.0,INSERT_VALUES); }
 55:     MatSetValue(A,i,i,2.0,INSERT_VALUES);
 56:   }
 57:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 58:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 59:   MatCreateVecs(A,&v,&w);
 60:   VecSet(v,1.0);

 62:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 63:                 Create the spectral transformation object
 64:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 66:   STCreate(PETSC_COMM_WORLD,&st);
 67:   mat[0] = A;
 68:   STSetOperators(st,1,mat);
 69:   STSetTransform(st,PETSC_TRUE);
 70:   STSetFromOptions(st);

 72:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 73:               Apply the transformed operator for several ST's
 74:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 76:   /* shift, sigma=0.0 */
 77:   STSetUp(st);
 78:   STGetType(st,&type);
 79:   PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type);
 80:   STApply(st,v,w);
 81:   VecView(w,NULL);

 83:   /* shift, sigma=0.1 */
 84:   sigma = 0.1;
 85:   STSetShift(st,sigma);
 86:   STGetShift(st,&sigma);
 87:   PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma));
 88:   STApply(st,v,w);
 89:   VecView(w,NULL);

 91:   /* sinvert, sigma=0.1 */
 92:   STPostSolve(st);   /* undo changes if inplace */
 93:   STSetType(st,STSINVERT);
 94:   STGetType(st,&type);
 95:   PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type);
 96:   STGetShift(st,&sigma);
 97:   PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma));
 98:   STApply(st,v,w);
 99:   VecView(w,NULL);

101:   /* sinvert, sigma=-0.5 */
102:   sigma = -0.5;
103:   STSetShift(st,sigma);
104:   STGetShift(st,&sigma);
105:   PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma));
106:   STApply(st,v,w);
107:   VecView(w,NULL);

109:   /* cayley, sigma=-0.5, tau=-0.5 (equal to sigma by default) */
110:   STPostSolve(st);   /* undo changes if inplace */
111:   STSetType(st,STCAYLEY);
112:   STSetUp(st);
113:   STGetType(st,&type);
114:   PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type);
115:   STGetShift(st,&sigma);
116:   STCayleyGetAntishift(st,&tau);
117:   PetscPrintf(PETSC_COMM_WORLD,"With shift=%g, antishift=%g\n",(double)PetscRealPart(sigma),(double)PetscRealPart(tau));
118:   STApply(st,v,w);
119:   VecView(w,NULL);

121:   /* cayley, sigma=1.1, tau=1.1 (still equal to sigma) */
122:   sigma = 1.1;
123:   STSetShift(st,sigma);
124:   STGetShift(st,&sigma);
125:   STCayleyGetAntishift(st,&tau);
126:   PetscPrintf(PETSC_COMM_WORLD,"With shift=%g, antishift=%g\n",(double)PetscRealPart(sigma),(double)PetscRealPart(tau));
127:   STApply(st,v,w);
128:   VecView(w,NULL);

130:   /* cayley, sigma=1.1, tau=-1.0 */
131:   tau = -1.0;
132:   STCayleySetAntishift(st,tau);
133:   STGetShift(st,&sigma);
134:   STCayleyGetAntishift(st,&tau);
135:   PetscPrintf(PETSC_COMM_WORLD,"With shift=%g, antishift=%g\n",(double)PetscRealPart(sigma),(double)PetscRealPart(tau));
136:   STApply(st,v,w);
137:   VecView(w,NULL);

139:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
140:                   Check inner product matrix in Cayley
141:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
142:   STGetBilinearForm(st,&B);
143:   MatMult(B,v,w);
144:   VecView(w,NULL);

146:   STDestroy(&st);
147:   MatDestroy(&A);
148:   MatDestroy(&B);
149:   VecDestroy(&v);
150:   VecDestroy(&w);
151:   SlepcFinalize();
152:   return ierr;
153: }