Actual source code: test3.c
slepc-3.18.3 2023-03-24
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;
25: SlepcInitialize(&argc,&argv,(char*)0,help);
26: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
27: 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: MatCreate(PETSC_COMM_WORLD,&A);
34: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
35: MatSetFromOptions(A);
36: MatSetUp(A);
38: MatCreate(PETSC_COMM_WORLD,&B);
39: MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n);
40: MatSetFromOptions(B);
41: MatSetUp(B);
43: MatGetOwnershipRange(A,&Istart,&Iend);
44: for (i=Istart;i<Iend;i++) {
45: MatSetValue(A,i,i,2.0,INSERT_VALUES);
46: if (i>0) {
47: MatSetValue(A,i,i-1,-1.0,INSERT_VALUES);
48: MatSetValue(B,i,i,(PetscScalar)i,INSERT_VALUES);
49: } else MatSetValue(B,i,i,-1.0,INSERT_VALUES);
50: if (i<n-1) MatSetValue(A,i,i+1,-1.0,INSERT_VALUES);
51: }
52: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
53: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
54: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
55: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
56: MatCreateVecs(A,&v,&w);
57: VecSet(v,1.0);
59: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
60: Create the spectral transformation object
61: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
63: STCreate(PETSC_COMM_WORLD,&st);
64: mat[0] = A;
65: mat[1] = B;
66: STSetMatrices(st,2,mat);
67: STSetTransform(st,PETSC_TRUE);
68: STSetFromOptions(st);
70: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
71: Apply the transformed operator for several ST's
72: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
74: /* shift, sigma=0.0 */
75: STSetUp(st);
76: STGetType(st,&type);
77: PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type);
78: STApply(st,v,w);
79: VecView(w,NULL);
81: /* shift, sigma=0.1 */
82: sigma = 0.1;
83: STSetShift(st,sigma);
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.1 */
90: STPostSolve(st); /* undo changes if inplace */
91: STSetType(st,STSINVERT);
92: STGetType(st,&type);
93: PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type);
94: STGetShift(st,&sigma);
95: PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma));
96: STApply(st,v,w);
97: VecView(w,NULL);
99: /* sinvert, sigma=-0.5 */
100: sigma = -0.5;
101: STSetShift(st,sigma);
102: STGetShift(st,&sigma);
103: PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma));
104: STApply(st,v,w);
105: VecView(w,NULL);
107: /* cayley, sigma=-0.5, tau=-0.5 (equal to sigma by default) */
108: STPostSolve(st); /* undo changes if inplace */
109: STSetType(st,STCAYLEY);
110: STSetUp(st);
111: STGetType(st,&type);
112: PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type);
113: STGetShift(st,&sigma);
114: STCayleyGetAntishift(st,&tau);
115: PetscPrintf(PETSC_COMM_WORLD,"With shift=%g, antishift=%g\n",(double)PetscRealPart(sigma),(double)PetscRealPart(tau));
116: STApply(st,v,w);
117: VecView(w,NULL);
119: /* cayley, sigma=1.1, tau=1.1 (still equal to sigma) */
120: sigma = 1.1;
121: STSetShift(st,sigma);
122: STGetShift(st,&sigma);
123: STCayleyGetAntishift(st,&tau);
124: PetscPrintf(PETSC_COMM_WORLD,"With shift=%g, antishift=%g\n",(double)PetscRealPart(sigma),(double)PetscRealPart(tau));
125: STApply(st,v,w);
126: VecView(w,NULL);
128: /* cayley, sigma=1.1, tau=-1.0 */
129: tau = -1.0;
130: STCayleySetAntishift(st,tau);
131: STGetShift(st,&sigma);
132: STCayleyGetAntishift(st,&tau);
133: PetscPrintf(PETSC_COMM_WORLD,"With shift=%g, antishift=%g\n",(double)PetscRealPart(sigma),(double)PetscRealPart(tau));
134: STApply(st,v,w);
135: VecView(w,NULL);
137: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
138: Check inner product matrix in Cayley
139: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
140: STGetBilinearForm(st,&M);
141: MatMult(M,v,w);
142: VecView(w,NULL);
144: STDestroy(&st);
145: MatDestroy(&A);
146: MatDestroy(&B);
147: MatDestroy(&M);
148: VecDestroy(&v);
149: VecDestroy(&w);
150: SlepcFinalize();
151: return 0;
152: }
154: /*TEST
156: test:
157: suffix: 1
158: args: -st_matmode {{copy inplace shell}}
159: output_file: output/test3_1.out
160: requires: !single
162: TEST*/