Actual source code: test8.c
slepc-3.17.1 2022-04-11
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 and split preconditioner.\n\n";
13: #include <slepcst.h>
15: int main(int argc,char **argv)
16: {
17: Mat A,B,Pa,Pb,Pmat,mat[2];
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,"\n1-D Laplacian plus diagonal, n=%" PetscInt_FMT "\n\n",n);
30: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
31: Compute the operator matrices (1-D Laplacian and diagonal)
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: MatGetOwnershipRange(A,&Istart,&Iend);
45: for (i=Istart;i<Iend;i++) {
46: MatSetValue(A,i,i,2.0,INSERT_VALUES);
47: if (i>0) {
48: MatSetValue(A,i,i-1,-1.0,INSERT_VALUES);
49: MatSetValue(B,i,i,(PetscScalar)i,INSERT_VALUES);
50: } else MatSetValue(B,i,i,-1.0,INSERT_VALUES);
51: if (i<n-1) MatSetValue(A,i,i+1,-1.0,INSERT_VALUES);
52: }
53: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
54: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
55: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
56: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
57: MatCreateVecs(A,&v,&w);
58: VecSet(v,1.0);
60: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
61: Compute the split preconditioner matrices (two diagonals)
62: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
64: MatCreate(PETSC_COMM_WORLD,&Pa);
65: MatSetSizes(Pa,PETSC_DECIDE,PETSC_DECIDE,n,n);
66: MatSetFromOptions(Pa);
67: MatSetUp(Pa);
69: MatCreate(PETSC_COMM_WORLD,&Pb);
70: MatSetSizes(Pb,PETSC_DECIDE,PETSC_DECIDE,n,n);
71: MatSetFromOptions(Pb);
72: MatSetUp(Pb);
74: MatGetOwnershipRange(Pa,&Istart,&Iend);
75: for (i=Istart;i<Iend;i++) {
76: MatSetValue(Pa,i,i,2.0,INSERT_VALUES);
77: if (i>0) MatSetValue(Pb,i,i,(PetscScalar)i,INSERT_VALUES);
78: else MatSetValue(Pb,i,i,-1.0,INSERT_VALUES);
79: }
80: MatAssemblyBegin(Pa,MAT_FINAL_ASSEMBLY);
81: MatAssemblyEnd(Pa,MAT_FINAL_ASSEMBLY);
82: MatAssemblyBegin(Pb,MAT_FINAL_ASSEMBLY);
83: MatAssemblyEnd(Pb,MAT_FINAL_ASSEMBLY);
85: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
86: Create the spectral transformation object
87: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
89: STCreate(PETSC_COMM_WORLD,&st);
90: mat[0] = A;
91: mat[1] = B;
92: STSetMatrices(st,2,mat);
93: mat[0] = Pa;
94: mat[1] = Pb;
95: STSetSplitPreconditioner(st,2,mat,SAME_NONZERO_PATTERN);
96: STSetTransform(st,PETSC_TRUE);
97: STSetFromOptions(st);
98: STCayleySetAntishift(st,-0.2); /* only relevant for cayley */
100: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
101: Form the preconditioner matrix and print it
102: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
104: STGetKSP(st,&ksp);
105: KSPGetPC(ksp,&pc);
106: STGetOperator(st,NULL);
107: PCGetOperators(pc,NULL,&Pmat);
108: MatView(Pmat,NULL);
110: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
111: Apply the operator
112: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
114: /* sigma=0.0 */
115: STSetUp(st);
116: STGetType(st,&type);
117: PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type);
118: STApply(st,v,w);
119: VecView(w,NULL);
121: /* sigma=0.1 */
122: sigma = 0.1;
123: STSetShift(st,sigma);
124: STGetShift(st,&sigma);
125: PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma));
126: STGetOperator(st,NULL);
127: PCGetOperators(pc,NULL,&Pmat);
128: MatView(Pmat,NULL);
129: STApply(st,v,w);
130: VecView(w,NULL);
132: STDestroy(&st);
133: MatDestroy(&A);
134: MatDestroy(&B);
135: MatDestroy(&Pa);
136: MatDestroy(&Pb);
137: VecDestroy(&v);
138: VecDestroy(&w);
139: SlepcFinalize();
140: return 0;
141: }
143: /*TEST
145: test:
146: suffix: 1
147: args: -st_type {{cayley shift sinvert}separate output}
148: requires: !single
150: TEST*/