Actual source code: test3.c

slepc-3.17.1 2022-04-11
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, 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 two matrices.\n\n";

 13: #include <slepcst.h>

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

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

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

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

 37:   MatCreate(PETSC_COMM_WORLD,&B);
 38:   MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n);
 39:   MatSetFromOptions(B);
 40:   MatSetUp(B);

 42:   MatGetOwnershipRange(A,&Istart,&Iend);
 43:   for (i=Istart;i<Iend;i++) {
 44:     MatSetValue(A,i,i,2.0,INSERT_VALUES);
 45:     if (i>0) {
 46:       MatSetValue(A,i,i-1,-1.0,INSERT_VALUES);
 47:       MatSetValue(B,i,i,(PetscScalar)i,INSERT_VALUES);
 48:     } else MatSetValue(B,i,i,-1.0,INSERT_VALUES);
 49:     if (i<n-1) MatSetValue(A,i,i+1,-1.0,INSERT_VALUES);
 50:   }
 51:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 52:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 53:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 54:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
 55:   MatCreateVecs(A,&v,&w);
 56:   VecSet(v,1.0);

 58:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 59:                 Create the spectral transformation object
 60:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 62:   STCreate(PETSC_COMM_WORLD,&st);
 63:   mat[0] = A;
 64:   mat[1] = B;
 65:   STSetMatrices(st,2,mat);
 66:   STSetTransform(st,PETSC_TRUE);
 67:   STSetFromOptions(st);

 69:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 70:               Apply the transformed operator for several ST's
 71:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

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

 88:   /* sinvert, sigma=0.1 */
 89:   STPostSolve(st);   /* undo changes if inplace */
 90:   STSetType(st,STSINVERT);
 91:   STGetType(st,&type);
 92:   PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type);
 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:   /* sinvert, sigma=-0.5 */
 99:   sigma = -0.5;
100:   STSetShift(st,sigma);
101:   STGetShift(st,&sigma);
102:   PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma));
103:   STApply(st,v,w);
104:   VecView(w,NULL);

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

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

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

136:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
137:                   Check inner product matrix in Cayley
138:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
139:   STGetBilinearForm(st,&M);
140:   MatMult(M,v,w);
141:   VecView(w,NULL);

143:   STDestroy(&st);
144:   MatDestroy(&A);
145:   MatDestroy(&B);
146:   MatDestroy(&M);
147:   VecDestroy(&v);
148:   VecDestroy(&w);
149:   SlepcFinalize();
150:   return 0;
151: }

153: /*TEST

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

161: TEST*/