Actual source code: test2.c

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

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */

 11: static char help[] = "Test ST with one matrix.\n\n";

 13: #include <slepcst.h>

 15: int main(int argc,char **argv)
 16: {
 17:   Mat            A,B,mat[1];
 18:   ST             st;
 19:   Vec            v,w;
 20:   STType         type;
 21:   PetscScalar    sigma,tau;
 22:   PetscInt       n=10,i,Istart,Iend;

 25:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 26:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 27:   PetscPrintf(PETSC_COMM_WORLD,"\n1-D Laplacian, n=%D\n\n",n);

 29:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 30:      Compute the operator matrix for the 1-D Laplacian
 31:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 33:   MatCreate(PETSC_COMM_WORLD,&A);
 34:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
 35:   MatSetFromOptions(A);
 36:   MatSetUp(A);

 38:   MatGetOwnershipRange(A,&Istart,&Iend);
 39:   for (i=Istart;i<Iend;i++) {
 40:     if (i>0) { MatSetValue(A,i,i-1,-1.0,INSERT_VALUES); }
 41:     if (i<n-1) { MatSetValue(A,i,i+1,-1.0,INSERT_VALUES); }
 42:     MatSetValue(A,i,i,2.0,INSERT_VALUES);
 43:   }
 44:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 45:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 46:   MatCreateVecs(A,&v,&w);
 47:   VecSet(v,1.0);

 49:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 50:                 Create the spectral transformation object
 51:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 53:   STCreate(PETSC_COMM_WORLD,&st);
 54:   mat[0] = A;
 55:   STSetMatrices(st,1,mat);
 56:   STSetTransform(st,PETSC_TRUE);
 57:   STSetFromOptions(st);

 59:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 60:               Apply the transformed operator for several ST's
 61:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 63:   /* shift, sigma=0.0 */
 64:   STSetUp(st);
 65:   STGetType(st,&type);
 66:   PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type);
 67:   STApply(st,v,w);
 68:   VecView(w,NULL);

 70:   /* shift, sigma=0.1 */
 71:   sigma = 0.1;
 72:   STSetShift(st,sigma);
 73:   STGetShift(st,&sigma);
 74:   PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma));
 75:   STApply(st,v,w);
 76:   VecView(w,NULL);
 77:   STApplyTranspose(st,v,w);
 78:   VecView(w,NULL);

 80:   /* sinvert, sigma=0.1 */
 81:   STPostSolve(st);   /* undo changes if inplace */
 82:   STSetType(st,STSINVERT);
 83:   STGetType(st,&type);
 84:   PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type);
 85:   STGetShift(st,&sigma);
 86:   PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma));
 87:   STApply(st,v,w);
 88:   VecView(w,NULL);

 90:   /* sinvert, sigma=-0.5 */
 91:   sigma = -0.5;
 92:   STSetShift(st,sigma);
 93:   STGetShift(st,&sigma);
 94:   PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma));
 95:   STApply(st,v,w);
 96:   VecView(w,NULL);

 98:   /* cayley, sigma=-0.5, tau=-0.5 (equal to sigma by default) */
 99:   STPostSolve(st);   /* undo changes if inplace */
100:   STSetType(st,STCAYLEY);
101:   STSetUp(st);
102:   STGetType(st,&type);
103:   PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type);
104:   STGetShift(st,&sigma);
105:   STCayleyGetAntishift(st,&tau);
106:   PetscPrintf(PETSC_COMM_WORLD,"With shift=%g, antishift=%g\n",(double)PetscRealPart(sigma),(double)PetscRealPart(tau));
107:   STApply(st,v,w);
108:   VecView(w,NULL);

110:   /* cayley, sigma=1.1, tau=1.1 (still equal to sigma) */
111:   sigma = 1.1;
112:   STSetShift(st,sigma);
113:   STGetShift(st,&sigma);
114:   STCayleyGetAntishift(st,&tau);
115:   PetscPrintf(PETSC_COMM_WORLD,"With shift=%g, antishift=%g\n",(double)PetscRealPart(sigma),(double)PetscRealPart(tau));
116:   STApply(st,v,w);
117:   VecView(w,NULL);

119:   /* cayley, sigma=1.1, tau=-1.0 */
120:   tau = -1.0;
121:   STCayleySetAntishift(st,tau);
122:   STGetShift(st,&sigma);
123:   STCayleyGetAntishift(st,&tau);
124:   PetscPrintf(PETSC_COMM_WORLD,"With shift=%g, antishift=%g\n",(double)PetscRealPart(sigma),(double)PetscRealPart(tau));
125:   STApply(st,v,w);
126:   VecView(w,NULL);

128:   /* check inner product matrix in Cayley */
129:   STGetBilinearForm(st,&B);
130:   MatMult(B,v,w);
131:   VecView(w,NULL);

133:   /* check again sinvert, sigma=0.1 */
134:   STPostSolve(st);   /* undo changes if inplace */
135:   STSetType(st,STSINVERT);
136:   STGetType(st,&type);
137:   PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type);
138:   sigma = 0.1;
139:   STSetShift(st,sigma);
140:   STGetShift(st,&sigma);
141:   PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma));
142:   STApply(st,v,w);
143:   VecView(w,NULL);

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

154: /*TEST

156:    test:
157:       suffix: 1
158:       args: -st_matmode {{copy inplace shell}}
159:       output_file: output/test2_1.out
160:       requires: !single

162: TEST*/