Actual source code: test2.c

slepc-3.20.1 2023-11-27
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:   PetscFunctionBeginUser;
 25:   PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
 26:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 27:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n1-D Laplacian, n=%" PetscInt_FMT "\n\n",n));

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

145:   PetscCall(STDestroy(&st));
146:   PetscCall(MatDestroy(&A));
147:   PetscCall(MatDestroy(&B));
148:   PetscCall(VecDestroy(&v));
149:   PetscCall(VecDestroy(&w));
150:   PetscCall(SlepcFinalize());
151:   return 0;
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*/