Actual source code: test4.c
slepc-3.17.0 2022-03-31
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 four matrices.\n\n";
13: #include <slepcst.h>
15: int main(int argc,char **argv)
16: {
17: Mat A,B,C,D,mat[4];
18: ST st;
19: KSP ksp;
20: Vec v,w;
21: STType type;
22: PetscScalar sigma;
23: 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,"\nTest ST with four matrices, n=%" PetscInt_FMT "\n\n",n);
28: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
29: Compute the operator matrices
30: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
32: MatCreate(PETSC_COMM_WORLD,&A);
33: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
34: MatSetFromOptions(A);
35: MatSetUp(A);
37: MatCreate(PETSC_COMM_WORLD,&B);
38: MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n);
39: MatSetFromOptions(B);
40: MatSetUp(B);
42: MatCreate(PETSC_COMM_WORLD,&C);
43: MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n);
44: MatSetFromOptions(C);
45: MatSetUp(C);
47: MatCreate(PETSC_COMM_WORLD,&D);
48: MatSetSizes(D,PETSC_DECIDE,PETSC_DECIDE,n,n);
49: MatSetFromOptions(D);
50: MatSetUp(D);
52: MatGetOwnershipRange(A,&Istart,&Iend);
53: for (i=Istart;i<Iend;i++) {
54: MatSetValue(A,i,i,2.0,INSERT_VALUES);
55: if (i>0) {
56: MatSetValue(A,i,i-1,-1.0,INSERT_VALUES);
57: MatSetValue(B,i,i,(PetscScalar)i,INSERT_VALUES);
58: } else MatSetValue(B,i,i,-1.0,INSERT_VALUES);
59: if (i<n-1) MatSetValue(A,i,i+1,-1.0,INSERT_VALUES);
60: MatSetValue(C,i,n-i-1,1.0,INSERT_VALUES);
61: MatSetValue(D,i,i,i*.1,INSERT_VALUES);
62: if (i==0) MatSetValue(D,0,n-1,1.0,INSERT_VALUES);
63: if (i==n-1) MatSetValue(D,n-1,0,1.0,INSERT_VALUES);
64: }
66: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
67: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
68: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
69: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
70: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
71: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
72: MatAssemblyBegin(D,MAT_FINAL_ASSEMBLY);
73: MatAssemblyEnd(D,MAT_FINAL_ASSEMBLY);
74: MatCreateVecs(A,&v,&w);
75: VecSet(v,1.0);
77: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
78: Create the spectral transformation object
79: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
81: STCreate(PETSC_COMM_WORLD,&st);
82: mat[0] = A;
83: mat[1] = B;
84: mat[2] = C;
85: mat[3] = D;
86: STSetMatrices(st,4,mat);
87: STGetKSP(st,&ksp);
88: KSPSetTolerances(ksp,100*PETSC_MACHINE_EPSILON,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
89: STSetFromOptions(st);
91: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
92: Apply the transformed operator for several ST's
93: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
95: /* shift, sigma=0.0 */
96: STSetUp(st);
97: STGetType(st,&type);
98: PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type);
99: for (i=0;i<4;i++) {
100: STMatMult(st,i,v,w);
101: PetscPrintf(PETSC_COMM_WORLD,"k= %" PetscInt_FMT "\n",i);
102: VecView(w,NULL);
103: }
104: STMatSolve(st,v,w);
105: PetscPrintf(PETSC_COMM_WORLD,"solve\n");
106: VecView(w,NULL);
108: /* shift, sigma=0.1 */
109: sigma = 0.1;
110: STSetShift(st,sigma);
111: STGetShift(st,&sigma);
112: PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma));
113: for (i=0;i<4;i++) {
114: STMatMult(st,i,v,w);
115: PetscPrintf(PETSC_COMM_WORLD,"k= %" PetscInt_FMT "\n",i);
116: VecView(w,NULL);
117: }
118: STMatSolve(st,v,w);
119: PetscPrintf(PETSC_COMM_WORLD,"solve\n");
120: VecView(w,NULL);
122: /* sinvert, sigma=0.1 */
123: STPostSolve(st);
124: STSetType(st,STSINVERT);
125: STGetType(st,&type);
126: PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type);
127: STGetShift(st,&sigma);
128: PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma));
129: for (i=0;i<4;i++) {
130: STMatMult(st,i,v,w);
131: PetscPrintf(PETSC_COMM_WORLD,"k= %" PetscInt_FMT "\n",i);
132: VecView(w,NULL);
133: }
134: STMatSolve(st,v,w);
135: PetscPrintf(PETSC_COMM_WORLD,"solve\n");
136: VecView(w,NULL);
138: /* sinvert, sigma=-0.5 */
139: sigma = -0.5;
140: STSetShift(st,sigma);
141: STGetShift(st,&sigma);
142: PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma));
143: for (i=0;i<4;i++) {
144: STMatMult(st,i,v,w);
145: PetscPrintf(PETSC_COMM_WORLD,"k= %" PetscInt_FMT "\n",i);
146: VecView(w,NULL);
147: }
148: STMatSolve(st,v,w);
149: PetscPrintf(PETSC_COMM_WORLD,"solve\n");
150: VecView(w,NULL);
152: STDestroy(&st);
153: MatDestroy(&A);
154: MatDestroy(&B);
155: MatDestroy(&C);
156: MatDestroy(&D);
157: VecDestroy(&v);
158: VecDestroy(&w);
159: SlepcFinalize();
160: return 0;
161: }
163: /*TEST
165: test:
166: suffix: 1
167: args: -st_transform -st_matmode {{copy shell}}
168: output_file: output/test4_1.out
169: requires: !single
171: test:
172: suffix: 2
173: args: -st_matmode {{copy shell}}
174: output_file: output/test4_2.out
176: TEST*/