Actual source code: test4.c
slepc-3.7.4 2017-05-17
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
8: SLEPc is free software: you can redistribute it and/or modify it under the
9: terms of version 3 of the GNU Lesser General Public License as published by
10: the Free Software Foundation.
12: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
13: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
14: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
15: more details.
17: You should have received a copy of the GNU Lesser General Public License
18: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
19: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
20: */
22: static char help[] = "Test ST with four matrices.\n\n";
24: #include <slepcst.h>
28: int main(int argc,char **argv)
29: {
30: Mat A,B,C,D,mat[4];
31: ST st;
32: KSP ksp;
33: Vec v,w;
34: STType type;
35: PetscScalar sigma;
36: PetscInt n=10,i,Istart,Iend;
39: SlepcInitialize(&argc,&argv,(char*)0,help);
40: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
41: PetscPrintf(PETSC_COMM_WORLD,"\n1-D Laplacian plus diagonal, n=%D\n\n",n);
42: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
43: Compute the operator matrix for the 1-D Laplacian
44: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
46: MatCreate(PETSC_COMM_WORLD,&A);
47: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
48: MatSetFromOptions(A);
49: MatSetUp(A);
51: MatCreate(PETSC_COMM_WORLD,&B);
52: MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n);
53: MatSetFromOptions(B);
54: MatSetUp(B);
56: MatCreate(PETSC_COMM_WORLD,&C);
57: MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n);
58: MatSetFromOptions(C);
59: MatSetUp(C);
61: MatCreate(PETSC_COMM_WORLD,&D);
62: MatSetSizes(D,PETSC_DECIDE,PETSC_DECIDE,n,n);
63: MatSetFromOptions(D);
64: MatSetUp(D);
66: MatGetOwnershipRange(A,&Istart,&Iend);
67: for (i=Istart;i<Iend;i++) {
68: MatSetValue(A,i,i,2.0,INSERT_VALUES);
69: if (i>0) {
70: MatSetValue(A,i,i-1,-1.0,INSERT_VALUES);
71: MatSetValue(B,i,i,(PetscScalar)i,INSERT_VALUES);
72: } else {
73: MatSetValue(B,i,i,-1.0,INSERT_VALUES);
74: }
75: if (i<n-1) {
76: MatSetValue(A,i,i+1,-1.0,INSERT_VALUES);
77: }
78: MatSetValue(C,i,n-i-1,1.0,INSERT_VALUES);
79: MatSetValue(D,i,i,i*.1,INSERT_VALUES);
80: if (i==0) {
81: MatSetValue(D,0,n-1,1.0,INSERT_VALUES);
82: }
83: if (i==n-1) {
84: MatSetValue(D,n-1,0,1.0,INSERT_VALUES);
85: }
86: }
88: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
89: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
90: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
91: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
92: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
93: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
94: MatAssemblyBegin(D,MAT_FINAL_ASSEMBLY);
95: MatAssemblyEnd(D,MAT_FINAL_ASSEMBLY);
96: MatCreateVecs(A,&v,&w);
97: VecSet(v,1.0);
99: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
100: Create the spectral transformation object
101: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
102: STCreate(PETSC_COMM_WORLD,&st);
103: mat[0] = A;
104: mat[1] = B;
105: mat[2] = C;
106: mat[3] = D;
107: STSetOperators(st,4,mat);
108: STGetKSP(st,&ksp);
109: KSPSetTolerances(ksp,100*PETSC_MACHINE_EPSILON,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
110: STSetFromOptions(st);
111: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
112: Apply the transformed operator for several ST's
113: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
114: /* shift, sigma=0.0 */
115: STSetUp(st);
116: STGetType(st,&type);
117: PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type);
118: for (i=0;i<4;i++) {
119: STMatMult(st,i,v,w);
120: PetscPrintf(PETSC_COMM_WORLD,"k= %D\n",i);
121: VecView(w,NULL);
122: }
123: STMatSolve(st,v,w);
124: PetscPrintf(PETSC_COMM_WORLD,"solve\n");
125: VecView(w,NULL);
127: /* shift, sigma=0.1 */
128: sigma = 0.1;
129: STSetShift(st,sigma);
130: STGetShift(st,&sigma);
131: PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma));
132: for (i=0;i<4;i++) {
133: STMatMult(st,i,v,w);
134: PetscPrintf(PETSC_COMM_WORLD,"k= %D\n",i);
135: VecView(w,NULL);
136: }
137: STMatSolve(st,v,w);
138: PetscPrintf(PETSC_COMM_WORLD,"solve\n");
139: VecView(w,NULL);
141: /* sinvert, sigma=0.1 */
142: STPostSolve(st);
143: STSetType(st,STSINVERT);
144: STGetType(st,&type);
145: PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type);
146: STGetShift(st,&sigma);
147: PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma));
148: for (i=0;i<4;i++) {
149: STMatMult(st,i,v,w);
150: PetscPrintf(PETSC_COMM_WORLD,"k= %D\n",i);
151: VecView(w,NULL);
152: }
153: STMatSolve(st,v,w);
154: PetscPrintf(PETSC_COMM_WORLD,"solve\n");
155: VecView(w,NULL);
157: /* sinvert, sigma=-0.5 */
158: sigma = -0.5;
159: STSetShift(st,sigma);
160: STGetShift(st,&sigma);
161: PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma));
162: for (i=0;i<4;i++) {
163: STMatMult(st,i,v,w);
164: PetscPrintf(PETSC_COMM_WORLD,"k= %D\n",i);
165: VecView(w,NULL);
166: }
167: STMatSolve(st,v,w);
168: PetscPrintf(PETSC_COMM_WORLD,"solve\n");
169: VecView(w,NULL);
170: STDestroy(&st);
171: MatDestroy(&A);
172: MatDestroy(&B);
173: MatDestroy(&C);
174: MatDestroy(&D);
175: VecDestroy(&v);
176: VecDestroy(&w);
177: SlepcFinalize();
178: return ierr;
179: }