Actual source code: test8.c
slepc-3.20.2 2024-03-15
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: PetscFunctionBeginUser;
27: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
28: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
29: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n1-D Laplacian plus diagonal, n=%" PetscInt_FMT "\n\n",n));
31: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
32: Compute the operator matrices (1-D Laplacian and diagonal)
33: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
35: PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
36: PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n));
37: PetscCall(MatSetFromOptions(A));
38: PetscCall(MatSetUp(A));
40: PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
41: PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n));
42: PetscCall(MatSetFromOptions(B));
43: PetscCall(MatSetUp(B));
45: PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
46: for (i=Istart;i<Iend;i++) {
47: PetscCall(MatSetValue(A,i,i,2.0,INSERT_VALUES));
48: if (i>0) {
49: PetscCall(MatSetValue(A,i,i-1,-1.0,INSERT_VALUES));
50: PetscCall(MatSetValue(B,i,i,(PetscScalar)i,INSERT_VALUES));
51: } else PetscCall(MatSetValue(B,i,i,-1.0,INSERT_VALUES));
52: if (i<n-1) PetscCall(MatSetValue(A,i,i+1,-1.0,INSERT_VALUES));
53: }
54: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
55: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
56: PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
57: PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
58: PetscCall(MatCreateVecs(A,&v,&w));
59: PetscCall(VecSet(v,1.0));
61: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
62: Compute the split preconditioner matrices (two diagonals)
63: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
65: PetscCall(MatCreate(PETSC_COMM_WORLD,&Pa));
66: PetscCall(MatSetSizes(Pa,PETSC_DECIDE,PETSC_DECIDE,n,n));
67: PetscCall(MatSetFromOptions(Pa));
68: PetscCall(MatSetUp(Pa));
70: PetscCall(MatCreate(PETSC_COMM_WORLD,&Pb));
71: PetscCall(MatSetSizes(Pb,PETSC_DECIDE,PETSC_DECIDE,n,n));
72: PetscCall(MatSetFromOptions(Pb));
73: PetscCall(MatSetUp(Pb));
75: PetscCall(MatGetOwnershipRange(Pa,&Istart,&Iend));
76: for (i=Istart;i<Iend;i++) {
77: PetscCall(MatSetValue(Pa,i,i,2.0,INSERT_VALUES));
78: if (i>0) PetscCall(MatSetValue(Pb,i,i,(PetscScalar)i,INSERT_VALUES));
79: else PetscCall(MatSetValue(Pb,i,i,-1.0,INSERT_VALUES));
80: }
81: PetscCall(MatAssemblyBegin(Pa,MAT_FINAL_ASSEMBLY));
82: PetscCall(MatAssemblyEnd(Pa,MAT_FINAL_ASSEMBLY));
83: PetscCall(MatAssemblyBegin(Pb,MAT_FINAL_ASSEMBLY));
84: PetscCall(MatAssemblyEnd(Pb,MAT_FINAL_ASSEMBLY));
86: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
87: Create the spectral transformation object
88: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
90: PetscCall(STCreate(PETSC_COMM_WORLD,&st));
91: mat[0] = A;
92: mat[1] = B;
93: PetscCall(STSetMatrices(st,2,mat));
94: mat[0] = Pa;
95: mat[1] = Pb;
96: PetscCall(STSetSplitPreconditioner(st,2,mat,SAME_NONZERO_PATTERN));
97: PetscCall(STSetTransform(st,PETSC_TRUE));
98: PetscCall(STSetFromOptions(st));
99: PetscCall(STCayleySetAntishift(st,-0.2)); /* only relevant for cayley */
101: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
102: Form the preconditioner matrix and print it
103: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
105: PetscCall(STGetKSP(st,&ksp));
106: PetscCall(KSPGetPC(ksp,&pc));
107: PetscCall(STGetOperator(st,NULL));
108: PetscCall(PCGetOperators(pc,NULL,&Pmat));
109: PetscCall(MatView(Pmat,NULL));
111: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
112: Apply the operator
113: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
115: /* sigma=0.0 */
116: PetscCall(STSetUp(st));
117: PetscCall(STGetType(st,&type));
118: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type));
119: PetscCall(STApply(st,v,w));
120: PetscCall(VecView(w,NULL));
122: /* sigma=0.1 */
123: sigma = 0.1;
124: PetscCall(STSetShift(st,sigma));
125: PetscCall(STGetShift(st,&sigma));
126: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma)));
127: PetscCall(STGetOperator(st,NULL));
128: PetscCall(PCGetOperators(pc,NULL,&Pmat));
129: PetscCall(MatView(Pmat,NULL));
130: PetscCall(STApply(st,v,w));
131: PetscCall(VecView(w,NULL));
133: PetscCall(STDestroy(&st));
134: PetscCall(MatDestroy(&A));
135: PetscCall(MatDestroy(&B));
136: PetscCall(MatDestroy(&Pa));
137: PetscCall(MatDestroy(&Pb));
138: PetscCall(VecDestroy(&v));
139: PetscCall(VecDestroy(&w));
140: PetscCall(SlepcFinalize());
141: return 0;
142: }
144: /*TEST
146: test:
147: suffix: 1
148: args: -st_type {{cayley shift sinvert}separate output}
149: requires: !single
151: TEST*/