Actual source code: test7.c
slepc-3.19.2 2023-09-05
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));
39: PetscCall(MatSetUp(A));
41: PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
42: for (i=Istart;i<Iend;i++) {
43: if (i>0) PetscCall(MatSetValue(A,i,i-1,-1.0,INSERT_VALUES));
44: if (i<n-1) PetscCall(MatSetValue(A,i,i+1,-1.0,INSERT_VALUES));
45: PetscCall(MatSetValue(A,i,i,2.0,INSERT_VALUES));
46: }
47: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
48: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
49: PetscCall(MatCreateVecs(A,&v,&w));
50: PetscCall(VecSet(v,1.0));
52: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
53: Compute the split preconditioner matrix (one diagonal)
54: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
56: PetscCall(MatCreate(PETSC_COMM_WORLD,&Pa));
57: PetscCall(MatSetSizes(Pa,PETSC_DECIDE,PETSC_DECIDE,n,n));
58: PetscCall(MatSetFromOptions(Pa));
59: PetscCall(MatSetUp(Pa));
61: PetscCall(MatGetOwnershipRange(Pa,&Istart,&Iend));
62: for (i=Istart;i<Iend;i++) PetscCall(MatSetValue(Pa,i,i,2.0,INSERT_VALUES));
63: PetscCall(MatAssemblyBegin(Pa,MAT_FINAL_ASSEMBLY));
64: PetscCall(MatAssemblyEnd(Pa,MAT_FINAL_ASSEMBLY));
66: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
67: Create the spectral transformation object
68: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
70: PetscCall(STCreate(PETSC_COMM_WORLD,&st));
71: mat[0] = A;
72: PetscCall(STSetMatrices(st,1,mat));
73: mat[0] = Pa;
74: PetscCall(STSetSplitPreconditioner(st,1,mat,SAME_NONZERO_PATTERN));
75: PetscCall(STSetTransform(st,PETSC_TRUE));
76: PetscCall(STSetFromOptions(st));
77: PetscCall(STCayleySetAntishift(st,-0.3)); /* only relevant for cayley */
79: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
80: Form the preconditioner matrix and print it
81: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
83: PetscCall(PetscObjectTypeCompareAny((PetscObject)st,&flg,STSINVERT,STCAYLEY,""));
84: if (flg) {
85: PetscCall(STGetKSP(st,&ksp));
86: PetscCall(KSPGetPC(ksp,&pc));
87: PetscCall(STGetOperator(st,NULL));
88: PetscCall(PCGetOperators(pc,NULL,&Pmat));
89: PetscCall(MatView(Pmat,NULL));
90: }
92: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
93: Apply the operator
94: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
96: /* sigma=0.0 */
97: PetscCall(STSetUp(st));
98: PetscCall(STGetType(st,&type));
99: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type));
100: PetscCall(STApply(st,v,w));
101: PetscCall(VecView(w,NULL));
103: /* sigma=0.1 */
104: sigma = 0.1;
105: PetscCall(STSetShift(st,sigma));
106: PetscCall(STGetShift(st,&sigma));
107: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma)));
108: if (flg) {
109: PetscCall(STGetOperator(st,NULL));
110: PetscCall(PCGetOperators(pc,NULL,&Pmat));
111: PetscCall(MatView(Pmat,NULL));
112: }
113: PetscCall(STApply(st,v,w));
114: PetscCall(VecView(w,NULL));
116: PetscCall(STDestroy(&st));
117: PetscCall(MatDestroy(&A));
118: PetscCall(MatDestroy(&Pa));
119: PetscCall(VecDestroy(&v));
120: PetscCall(VecDestroy(&w));
121: PetscCall(SlepcFinalize());
122: return 0;
123: }
125: /*TEST
127: test:
128: suffix: 1
129: args: -st_type {{cayley shift sinvert}separate output}
130: requires: !single
132: TEST*/