Actual source code: test2.c

slepc-3.21.1 2024-04-26
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,val;
 22:   PetscInt       n=10,i,Istart,Iend;
 23:   PetscBool      test_compl=PETSC_FALSE;

 25:   PetscFunctionBeginUser;
 26:   PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
 27:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 28:   PetscCall(PetscOptionsGetBool(NULL,NULL,"-complex",&test_compl,NULL));
 29:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n1-D Laplacian, n=%" PetscInt_FMT "\n\n",n));

 31:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 32:      Compute the operator matrix for the 1-D Laplacian
 33:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

 39:   PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
 40: #if defined(PETSC_USE_COMPLEX)
 41:   val = test_compl? PetscCMPLX(-1.0,0.4): -1.0;
 42: #else
 43:   val = -1.0;
 44: #endif
 45:   for (i=Istart;i<Iend;i++) {
 46:     if (i>0) PetscCall(MatSetValue(A,i,i-1,val,INSERT_VALUES));
 47:     if (i<n-1) PetscCall(MatSetValue(A,i,i+1,PetscConj(val),INSERT_VALUES));
 48:     PetscCall(MatSetValue(A,i,i,2.0,INSERT_VALUES));
 49:   }
 50:   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
 51:   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
 52:   PetscCall(MatCreateVecs(A,&v,&w));
 53: #if defined(PETSC_USE_COMPLEX)
 54:   val = test_compl? PetscCMPLX(1.0,-0.01): 1.0;
 55: #else
 56:   val = 1.0;
 57: #endif
 58:   PetscCall(VecSet(v,val));

 60:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 61:                 Create the spectral transformation object
 62:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

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

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

 81:   /* shift, sigma=0.1 */
 82:   sigma = 0.1;
 83:   PetscCall(STSetShift(st,sigma));
 84:   PetscCall(STGetShift(st,&sigma));
 85:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma)));
 86:   PetscCall(STApply(st,v,w));
 87:   PetscCall(VecView(w,NULL));
 88:   PetscCall(STApplyTranspose(st,v,w));
 89:   PetscCall(VecView(w,NULL));
 90:   if (test_compl) {
 91:     PetscCall(STApplyHermitianTranspose(st,v,w));
 92:     PetscCall(VecView(w,NULL));
 93:   }

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

105:   /* sinvert, sigma=-0.5 */
106:   sigma = -0.5;
107:   PetscCall(STSetShift(st,sigma));
108:   PetscCall(STGetShift(st,&sigma));
109:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma)));
110:   PetscCall(STApply(st,v,w));
111:   PetscCall(VecView(w,NULL));
112:   if (test_compl) {
113:     PetscCall(STApplyHermitianTranspose(st,v,w));
114:     PetscCall(VecView(w,NULL));
115:   }

117:   /* cayley, sigma=-0.5, tau=-0.5 (equal to sigma by default) */
118:   PetscCall(STPostSolve(st));   /* undo changes if inplace */
119:   PetscCall(STSetType(st,STCAYLEY));
120:   PetscCall(STSetUp(st));
121:   PetscCall(STGetType(st,&type));
122:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type));
123:   PetscCall(STGetShift(st,&sigma));
124:   PetscCall(STCayleyGetAntishift(st,&tau));
125:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"With shift=%g, antishift=%g\n",(double)PetscRealPart(sigma),(double)PetscRealPart(tau)));
126:   PetscCall(STApply(st,v,w));
127:   PetscCall(VecView(w,NULL));

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

138:   /* cayley, sigma=1.1, tau=-1.0 */
139:   tau = -1.0;
140:   PetscCall(STCayleySetAntishift(st,tau));
141:   PetscCall(STGetShift(st,&sigma));
142:   PetscCall(STCayleyGetAntishift(st,&tau));
143:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"With shift=%g, antishift=%g\n",(double)PetscRealPart(sigma),(double)PetscRealPart(tau)));
144:   PetscCall(STApply(st,v,w));
145:   PetscCall(VecView(w,NULL));

147:   /* check inner product matrix in Cayley */
148:   PetscCall(STGetBilinearForm(st,&B));
149:   PetscCall(MatMult(B,v,w));
150:   PetscCall(VecView(w,NULL));

152:   /* check again sinvert, sigma=0.1 */
153:   PetscCall(STPostSolve(st));   /* undo changes if inplace */
154:   PetscCall(STSetType(st,STSINVERT));
155:   PetscCall(STGetType(st,&type));
156:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type));
157:   sigma = 0.1;
158:   PetscCall(STSetShift(st,sigma));
159:   PetscCall(STGetShift(st,&sigma));
160:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma)));
161:   PetscCall(STApply(st,v,w));
162:   PetscCall(VecView(w,NULL));

164:   PetscCall(STDestroy(&st));
165:   PetscCall(MatDestroy(&A));
166:   PetscCall(MatDestroy(&B));
167:   PetscCall(VecDestroy(&v));
168:   PetscCall(VecDestroy(&w));
169:   PetscCall(SlepcFinalize());
170:   return 0;
171: }

173: /*TEST

175:    test:
176:       suffix: 1
177:       args: -st_matmode {{copy inplace shell}}
178:       requires: !single

180:    test:
181:       suffix: 2
182:       args: -complex -st_matmode {{copy inplace shell}}
183:       requires: complex !single

185: TEST*/