Actual source code: test9.c
slepc-3.18.1 2022-11-02
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;
27: SlepcInitialize(&argc,&argv,(char*)0,help);
28: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
29: PetscPrintf(PETSC_COMM_WORLD,"\nTest ST with four matrices, n=%" PetscInt_FMT "\n\n",n);
30: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
31: Compute the operator matrices
32: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
34: MatCreate(PETSC_COMM_WORLD,&A);
35: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
36: MatSetFromOptions(A);
37: MatSetUp(A);
39: MatCreate(PETSC_COMM_WORLD,&B);
40: MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n);
41: MatSetFromOptions(B);
42: MatSetUp(B);
44: MatCreate(PETSC_COMM_WORLD,&C);
45: MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n);
46: MatSetFromOptions(C);
47: MatSetUp(C);
49: MatCreate(PETSC_COMM_WORLD,&D);
50: MatSetSizes(D,PETSC_DECIDE,PETSC_DECIDE,n,n);
51: MatSetFromOptions(D);
52: MatSetUp(D);
54: MatGetOwnershipRange(A,&Istart,&Iend);
55: for (i=Istart;i<Iend;i++) {
56: MatSetValue(A,i,i,2.0,INSERT_VALUES);
57: if (i>0) {
58: MatSetValue(A,i,i-1,-1.0,INSERT_VALUES);
59: MatSetValue(B,i,i,(PetscScalar)i,INSERT_VALUES);
60: } else MatSetValue(B,i,i,-1.0,INSERT_VALUES);
61: if (i<n-1) MatSetValue(A,i,i+1,-1.0,INSERT_VALUES);
62: MatSetValue(C,i,n-i-1,1.0,INSERT_VALUES);
63: MatSetValue(D,i,i,i*.1,INSERT_VALUES);
64: if (i==0) MatSetValue(D,0,n-1,1.0,INSERT_VALUES);
65: if (i==n-1) MatSetValue(D,n-1,0,1.0,INSERT_VALUES);
66: }
68: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
69: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
70: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
71: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
72: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
73: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
74: MatAssemblyBegin(D,MAT_FINAL_ASSEMBLY);
75: MatAssemblyEnd(D,MAT_FINAL_ASSEMBLY);
76: MatCreateVecs(A,&v,&w);
77: VecSet(v,1.0);
79: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
80: Compute the split preconditioner matrices (four diagonals)
81: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
83: MatCreate(PETSC_COMM_WORLD,&Pa);
84: MatSetSizes(Pa,PETSC_DECIDE,PETSC_DECIDE,n,n);
85: MatSetFromOptions(Pa);
86: MatSetUp(Pa);
88: MatCreate(PETSC_COMM_WORLD,&Pb);
89: MatSetSizes(Pb,PETSC_DECIDE,PETSC_DECIDE,n,n);
90: MatSetFromOptions(Pb);
91: MatSetUp(Pb);
93: MatCreate(PETSC_COMM_WORLD,&Pc);
94: MatSetSizes(Pc,PETSC_DECIDE,PETSC_DECIDE,n,n);
95: MatSetFromOptions(Pc);
96: MatSetUp(Pc);
98: MatCreate(PETSC_COMM_WORLD,&Pd);
99: MatSetSizes(Pd,PETSC_DECIDE,PETSC_DECIDE,n,n);
100: MatSetFromOptions(Pd);
101: MatSetUp(Pd);
103: MatGetOwnershipRange(Pa,&Istart,&Iend);
104: for (i=Istart;i<Iend;i++) {
105: MatSetValue(Pa,i,i,2.0,INSERT_VALUES);
106: if (i>0) MatSetValue(Pb,i,i,(PetscScalar)i,INSERT_VALUES);
107: else MatSetValue(Pb,i,i,-1.0,INSERT_VALUES);
108: MatSetValue(Pd,i,i,i*.1,INSERT_VALUES);
109: }
111: MatAssemblyBegin(Pa,MAT_FINAL_ASSEMBLY);
112: MatAssemblyEnd(Pa,MAT_FINAL_ASSEMBLY);
113: MatAssemblyBegin(Pb,MAT_FINAL_ASSEMBLY);
114: MatAssemblyEnd(Pb,MAT_FINAL_ASSEMBLY);
115: MatAssemblyBegin(Pc,MAT_FINAL_ASSEMBLY);
116: MatAssemblyEnd(Pc,MAT_FINAL_ASSEMBLY);
117: MatAssemblyBegin(Pd,MAT_FINAL_ASSEMBLY);
118: MatAssemblyEnd(Pd,MAT_FINAL_ASSEMBLY);
120: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
121: Create the spectral transformation object
122: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
124: STCreate(PETSC_COMM_WORLD,&st);
125: mat[0] = A;
126: mat[1] = B;
127: mat[2] = C;
128: mat[3] = D;
129: STSetMatrices(st,4,mat);
130: mat[0] = Pa;
131: mat[1] = Pb;
132: mat[2] = Pc;
133: mat[3] = Pd;
134: STSetSplitPreconditioner(st,4,mat,SUBSET_NONZERO_PATTERN);
135: STGetKSP(st,&ksp);
136: KSPSetTolerances(ksp,100*PETSC_MACHINE_EPSILON,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
137: STSetTransform(st,PETSC_TRUE);
138: STSetFromOptions(st);
139: STGetKSP(st,&ksp);
140: KSPGetPC(ksp,&pc);
142: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
143: Apply the operator
144: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
146: /* sigma=0.0 */
147: STSetUp(st);
148: STGetType(st,&type);
149: PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type);
150: PCGetOperators(pc,NULL,&Pmat);
151: MatView(Pmat,NULL);
152: STMatSolve(st,v,w);
153: VecView(w,NULL);
155: /* sigma=0.1 */
156: sigma = 0.1;
157: STSetShift(st,sigma);
158: STGetShift(st,&sigma);
159: PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma));
160: PCGetOperators(pc,NULL,&Pmat);
161: MatView(Pmat,NULL);
162: STMatSolve(st,v,w);
163: VecView(w,NULL);
165: STDestroy(&st);
166: MatDestroy(&A);
167: MatDestroy(&B);
168: MatDestroy(&C);
169: MatDestroy(&D);
170: MatDestroy(&Pa);
171: MatDestroy(&Pb);
172: MatDestroy(&Pc);
173: MatDestroy(&Pd);
174: VecDestroy(&v);
175: VecDestroy(&w);
176: SlepcFinalize();
177: return 0;
178: }
180: /*TEST
182: test:
183: suffix: 1
184: args: -st_type {{shift sinvert}separate output} -st_pc_type jacobi
185: requires: !single
187: TEST*/