Actual source code: test9.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 and split preconditioner.\n\n";
13: #include <slepcst.h>
15: int main(int argc,char **argv)
16: {
17: Mat A,B,C,D,Pa,Pb,Pc,Pd,Pmat,mat[4];
18: ST st;
19: KSP ksp;
20: PC pc;
21: Vec v,w;
22: STType type;
23: PetscScalar sigma;
24: PetscInt n=10,i,Istart,Iend;
26: SlepcInitialize(&argc,&argv,(char*)0,help);
27: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
28: PetscPrintf(PETSC_COMM_WORLD,"\nTest ST with four matrices, n=%" PetscInt_FMT "\n\n",n);
29: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
30: Compute the operator matrices
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: MatCreate(PETSC_COMM_WORLD,&C);
44: MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n);
45: MatSetFromOptions(C);
46: MatSetUp(C);
48: MatCreate(PETSC_COMM_WORLD,&D);
49: MatSetSizes(D,PETSC_DECIDE,PETSC_DECIDE,n,n);
50: MatSetFromOptions(D);
51: MatSetUp(D);
53: MatGetOwnershipRange(A,&Istart,&Iend);
54: for (i=Istart;i<Iend;i++) {
55: MatSetValue(A,i,i,2.0,INSERT_VALUES);
56: if (i>0) {
57: MatSetValue(A,i,i-1,-1.0,INSERT_VALUES);
58: MatSetValue(B,i,i,(PetscScalar)i,INSERT_VALUES);
59: } else MatSetValue(B,i,i,-1.0,INSERT_VALUES);
60: if (i<n-1) MatSetValue(A,i,i+1,-1.0,INSERT_VALUES);
61: MatSetValue(C,i,n-i-1,1.0,INSERT_VALUES);
62: MatSetValue(D,i,i,i*.1,INSERT_VALUES);
63: if (i==0) MatSetValue(D,0,n-1,1.0,INSERT_VALUES);
64: if (i==n-1) MatSetValue(D,n-1,0,1.0,INSERT_VALUES);
65: }
67: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
68: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
69: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
70: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
71: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
72: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
73: MatAssemblyBegin(D,MAT_FINAL_ASSEMBLY);
74: MatAssemblyEnd(D,MAT_FINAL_ASSEMBLY);
75: MatCreateVecs(A,&v,&w);
76: VecSet(v,1.0);
78: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
79: Compute the split preconditioner matrices (four diagonals)
80: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
82: MatCreate(PETSC_COMM_WORLD,&Pa);
83: MatSetSizes(Pa,PETSC_DECIDE,PETSC_DECIDE,n,n);
84: MatSetFromOptions(Pa);
85: MatSetUp(Pa);
87: MatCreate(PETSC_COMM_WORLD,&Pb);
88: MatSetSizes(Pb,PETSC_DECIDE,PETSC_DECIDE,n,n);
89: MatSetFromOptions(Pb);
90: MatSetUp(Pb);
92: MatCreate(PETSC_COMM_WORLD,&Pc);
93: MatSetSizes(Pc,PETSC_DECIDE,PETSC_DECIDE,n,n);
94: MatSetFromOptions(Pc);
95: MatSetUp(Pc);
97: MatCreate(PETSC_COMM_WORLD,&Pd);
98: MatSetSizes(Pd,PETSC_DECIDE,PETSC_DECIDE,n,n);
99: MatSetFromOptions(Pd);
100: MatSetUp(Pd);
102: MatGetOwnershipRange(Pa,&Istart,&Iend);
103: for (i=Istart;i<Iend;i++) {
104: MatSetValue(Pa,i,i,2.0,INSERT_VALUES);
105: if (i>0) MatSetValue(Pb,i,i,(PetscScalar)i,INSERT_VALUES);
106: else MatSetValue(Pb,i,i,-1.0,INSERT_VALUES);
107: MatSetValue(Pd,i,i,i*.1,INSERT_VALUES);
108: }
110: MatAssemblyBegin(Pa,MAT_FINAL_ASSEMBLY);
111: MatAssemblyEnd(Pa,MAT_FINAL_ASSEMBLY);
112: MatAssemblyBegin(Pb,MAT_FINAL_ASSEMBLY);
113: MatAssemblyEnd(Pb,MAT_FINAL_ASSEMBLY);
114: MatAssemblyBegin(Pc,MAT_FINAL_ASSEMBLY);
115: MatAssemblyEnd(Pc,MAT_FINAL_ASSEMBLY);
116: MatAssemblyBegin(Pd,MAT_FINAL_ASSEMBLY);
117: MatAssemblyEnd(Pd,MAT_FINAL_ASSEMBLY);
119: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
120: Create the spectral transformation object
121: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
123: STCreate(PETSC_COMM_WORLD,&st);
124: mat[0] = A;
125: mat[1] = B;
126: mat[2] = C;
127: mat[3] = D;
128: STSetMatrices(st,4,mat);
129: mat[0] = Pa;
130: mat[1] = Pb;
131: mat[2] = Pc;
132: mat[3] = Pd;
133: STSetSplitPreconditioner(st,4,mat,SUBSET_NONZERO_PATTERN);
134: STGetKSP(st,&ksp);
135: KSPSetTolerances(ksp,100*PETSC_MACHINE_EPSILON,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
136: STSetTransform(st,PETSC_TRUE);
137: STSetFromOptions(st);
138: STGetKSP(st,&ksp);
139: KSPGetPC(ksp,&pc);
141: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
142: Apply the operator
143: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
145: /* sigma=0.0 */
146: STSetUp(st);
147: STGetType(st,&type);
148: PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type);
149: PCGetOperators(pc,NULL,&Pmat);
150: MatView(Pmat,NULL);
151: STMatSolve(st,v,w);
152: VecView(w,NULL);
154: /* sigma=0.1 */
155: sigma = 0.1;
156: STSetShift(st,sigma);
157: STGetShift(st,&sigma);
158: PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma));
159: PCGetOperators(pc,NULL,&Pmat);
160: MatView(Pmat,NULL);
161: STMatSolve(st,v,w);
162: VecView(w,NULL);
164: STDestroy(&st);
165: MatDestroy(&A);
166: MatDestroy(&B);
167: MatDestroy(&C);
168: MatDestroy(&D);
169: MatDestroy(&Pa);
170: MatDestroy(&Pb);
171: MatDestroy(&Pc);
172: MatDestroy(&Pd);
173: VecDestroy(&v);
174: VecDestroy(&w);
175: SlepcFinalize();
176: return 0;
177: }
179: /*TEST
181: test:
182: suffix: 1
183: args: -st_type {{shift sinvert}separate output} -st_pc_type jacobi
184: requires: !single
186: TEST*/