Actual source code: test3.c
slepc-3.21.1 2024-04-26
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: 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 plus diagonal, 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));
37: PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
38: PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n));
39: PetscCall(MatSetFromOptions(B));
41: PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
42: for (i=Istart;i<Iend;i++) {
43: PetscCall(MatSetValue(A,i,i,2.0,INSERT_VALUES));
44: if (i>0) {
45: PetscCall(MatSetValue(A,i,i-1,-1.0,INSERT_VALUES));
46: PetscCall(MatSetValue(B,i,i,(PetscScalar)i,INSERT_VALUES));
47: } else PetscCall(MatSetValue(B,i,i,-1.0,INSERT_VALUES));
48: if (i<n-1) PetscCall(MatSetValue(A,i,i+1,-1.0,INSERT_VALUES));
49: }
50: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
51: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
52: PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
53: PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
54: PetscCall(MatCreateVecs(A,&v,&w));
55: PetscCall(VecSet(v,1.0));
57: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
58: Create the spectral transformation object
59: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
61: PetscCall(STCreate(PETSC_COMM_WORLD,&st));
62: mat[0] = A;
63: mat[1] = B;
64: PetscCall(STSetMatrices(st,2,mat));
65: PetscCall(STSetTransform(st,PETSC_TRUE));
66: PetscCall(STSetFromOptions(st));
68: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
69: Apply the transformed operator for several ST's
70: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
72: /* shift, sigma=0.0 */
73: PetscCall(STSetUp(st));
74: PetscCall(STGetType(st,&type));
75: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type));
76: PetscCall(STApply(st,v,w));
77: PetscCall(VecView(w,NULL));
79: /* shift, sigma=0.1 */
80: sigma = 0.1;
81: PetscCall(STSetShift(st,sigma));
82: PetscCall(STGetShift(st,&sigma));
83: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma)));
84: PetscCall(STApply(st,v,w));
85: PetscCall(VecView(w,NULL));
87: /* sinvert, sigma=0.1 */
88: PetscCall(STPostSolve(st)); /* undo changes if inplace */
89: PetscCall(STSetType(st,STSINVERT));
90: PetscCall(STGetType(st,&type));
91: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type));
92: PetscCall(STGetShift(st,&sigma));
93: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma)));
94: PetscCall(STApply(st,v,w));
95: PetscCall(VecView(w,NULL));
97: /* sinvert, sigma=-0.5 */
98: sigma = -0.5;
99: PetscCall(STSetShift(st,sigma));
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: /* cayley, sigma=-0.5, tau=-0.5 (equal to sigma by default) */
106: PetscCall(STPostSolve(st)); /* undo changes if inplace */
107: PetscCall(STSetType(st,STCAYLEY));
108: PetscCall(STSetUp(st));
109: PetscCall(STGetType(st,&type));
110: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type));
111: PetscCall(STGetShift(st,&sigma));
112: PetscCall(STCayleyGetAntishift(st,&tau));
113: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"With shift=%g, antishift=%g\n",(double)PetscRealPart(sigma),(double)PetscRealPart(tau)));
114: PetscCall(STApply(st,v,w));
115: PetscCall(VecView(w,NULL));
117: /* cayley, sigma=1.1, tau=1.1 (still equal to sigma) */
118: sigma = 1.1;
119: PetscCall(STSetShift(st,sigma));
120: PetscCall(STGetShift(st,&sigma));
121: PetscCall(STCayleyGetAntishift(st,&tau));
122: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"With shift=%g, antishift=%g\n",(double)PetscRealPart(sigma),(double)PetscRealPart(tau)));
123: PetscCall(STApply(st,v,w));
124: PetscCall(VecView(w,NULL));
126: /* cayley, sigma=1.1, tau=-1.0 */
127: tau = -1.0;
128: PetscCall(STCayleySetAntishift(st,tau));
129: PetscCall(STGetShift(st,&sigma));
130: PetscCall(STCayleyGetAntishift(st,&tau));
131: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"With shift=%g, antishift=%g\n",(double)PetscRealPart(sigma),(double)PetscRealPart(tau)));
132: PetscCall(STApply(st,v,w));
133: PetscCall(VecView(w,NULL));
135: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
136: Check inner product matrix in Cayley
137: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
138: PetscCall(STGetBilinearForm(st,&M));
139: PetscCall(MatMult(M,v,w));
140: PetscCall(VecView(w,NULL));
142: PetscCall(STDestroy(&st));
143: PetscCall(MatDestroy(&A));
144: PetscCall(MatDestroy(&B));
145: PetscCall(MatDestroy(&M));
146: PetscCall(VecDestroy(&v));
147: PetscCall(VecDestroy(&w));
148: PetscCall(SlepcFinalize());
149: return 0;
150: }
152: /*TEST
154: test:
155: suffix: 1
156: args: -st_matmode {{copy inplace shell}}
157: requires: !single
159: TEST*/