Actual source code: test7.c
slepc-3.21.2 2024-09-25
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 one matrix and split preconditioner.\n\n";
13: #include <slepcst.h>
15: int main(int argc,char **argv)
16: {
17: Mat A,Pa,Pmat,mat[1];
18: ST st;
19: KSP ksp;
20: PC pc;
21: Vec v,w;
22: STType type;
23: PetscBool flg;
24: PetscScalar sigma;
25: PetscInt n=10,i,Istart,Iend;
27: PetscFunctionBeginUser;
28: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
29: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
30: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n1-D Laplacian, n=%" PetscInt_FMT "\n\n",n));
32: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
33: Compute the operator matrix for the 1-D Laplacian
34: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
36: PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
37: PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n));
38: PetscCall(MatSetFromOptions(A));
40: PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
41: for (i=Istart;i<Iend;i++) {
42: if (i>0) PetscCall(MatSetValue(A,i,i-1,-1.0,INSERT_VALUES));
43: if (i<n-1) PetscCall(MatSetValue(A,i,i+1,-1.0,INSERT_VALUES));
44: PetscCall(MatSetValue(A,i,i,2.0,INSERT_VALUES));
45: }
46: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
47: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
48: PetscCall(MatCreateVecs(A,&v,&w));
49: PetscCall(VecSet(v,1.0));
51: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
52: Compute the split preconditioner matrix (one diagonal)
53: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
55: PetscCall(MatCreate(PETSC_COMM_WORLD,&Pa));
56: PetscCall(MatSetSizes(Pa,PETSC_DECIDE,PETSC_DECIDE,n,n));
57: PetscCall(MatSetFromOptions(Pa));
59: PetscCall(MatGetOwnershipRange(Pa,&Istart,&Iend));
60: for (i=Istart;i<Iend;i++) PetscCall(MatSetValue(Pa,i,i,2.0,INSERT_VALUES));
61: PetscCall(MatAssemblyBegin(Pa,MAT_FINAL_ASSEMBLY));
62: PetscCall(MatAssemblyEnd(Pa,MAT_FINAL_ASSEMBLY));
64: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
65: Create the spectral transformation object
66: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
68: PetscCall(STCreate(PETSC_COMM_WORLD,&st));
69: mat[0] = A;
70: PetscCall(STSetMatrices(st,1,mat));
71: mat[0] = Pa;
72: PetscCall(STSetSplitPreconditioner(st,1,mat,SAME_NONZERO_PATTERN));
73: PetscCall(STSetTransform(st,PETSC_TRUE));
74: PetscCall(STSetFromOptions(st));
75: PetscCall(STCayleySetAntishift(st,-0.3)); /* only relevant for cayley */
77: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
78: Form the preconditioner matrix and print it
79: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
81: PetscCall(PetscObjectTypeCompareAny((PetscObject)st,&flg,STSINVERT,STCAYLEY,""));
82: if (flg) {
83: PetscCall(STGetKSP(st,&ksp));
84: PetscCall(KSPGetPC(ksp,&pc));
85: PetscCall(STGetOperator(st,NULL));
86: PetscCall(PCGetOperators(pc,NULL,&Pmat));
87: PetscCall(MatView(Pmat,NULL));
88: }
90: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
91: Apply the operator
92: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
94: /* sigma=0.0 */
95: PetscCall(STSetUp(st));
96: PetscCall(STGetType(st,&type));
97: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type));
98: PetscCall(STApply(st,v,w));
99: PetscCall(VecView(w,NULL));
101: /* sigma=0.1 */
102: sigma = 0.1;
103: PetscCall(STSetShift(st,sigma));
104: PetscCall(STGetShift(st,&sigma));
105: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma)));
106: if (flg) {
107: PetscCall(STGetOperator(st,NULL));
108: PetscCall(PCGetOperators(pc,NULL,&Pmat));
109: PetscCall(MatView(Pmat,NULL));
110: }
111: PetscCall(STApply(st,v,w));
112: PetscCall(VecView(w,NULL));
114: PetscCall(STDestroy(&st));
115: PetscCall(MatDestroy(&A));
116: PetscCall(MatDestroy(&Pa));
117: PetscCall(VecDestroy(&v));
118: PetscCall(VecDestroy(&w));
119: PetscCall(SlepcFinalize());
120: return 0;
121: }
123: /*TEST
125: test:
126: suffix: 1
127: args: -st_type {{cayley shift sinvert}separate output}
128: requires: !single
130: TEST*/