Actual source code: test2.c

slepc-3.17.2 2022-08-09
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 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;

 24:   SlepcInitialize(&argc,&argv,(char*)0,help);
 25:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 26:   PetscPrintf(PETSC_COMM_WORLD,"\n1-D Laplacian, 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:   MatGetOwnershipRange(A,&Istart,&Iend);
 38:   for (i=Istart;i<Iend;i++) {
 39:     if (i>0) MatSetValue(A,i,i-1,-1.0,INSERT_VALUES);
 40:     if (i<n-1) MatSetValue(A,i,i+1,-1.0,INSERT_VALUES);
 41:     MatSetValue(A,i,i,2.0,INSERT_VALUES);
 42:   }
 43:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 44:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 45:   MatCreateVecs(A,&v,&w);
 46:   VecSet(v,1.0);

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

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

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

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

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

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

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

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

109:   /* cayley, sigma=1.1, tau=1.1 (still equal to sigma) */
110:   sigma = 1.1;
111:   STSetShift(st,sigma);
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.0 */
119:   tau = -1.0;
120:   STCayleySetAntishift(st,tau);
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:   /* check inner product matrix in Cayley */
128:   STGetBilinearForm(st,&B);
129:   MatMult(B,v,w);
130:   VecView(w,NULL);

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

144:   STDestroy(&st);
145:   MatDestroy(&A);
146:   MatDestroy(&B);
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/test2_1.out
159:       requires: !single

161: TEST*/