Line data Source code
1 : /*
2 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3 : SLEPc - Scalable Library for Eigenvalue Problem Computations
4 : Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
5 :
6 : This file is part of SLEPc.
7 : SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9 : */
10 :
11 : static char help[] = "Test ST with one matrix and split preconditioner.\n\n";
12 :
13 : #include <slepcst.h>
14 :
15 3 : int main(int argc,char **argv)
16 : {
17 3 : Mat A,Pa,Pmat,mat[1];
18 3 : ST st;
19 3 : KSP ksp;
20 3 : PC pc;
21 3 : Vec v,w;
22 3 : STType type;
23 3 : PetscBool flg;
24 3 : PetscScalar sigma;
25 3 : PetscInt n=10,i,Istart,Iend;
26 :
27 3 : PetscFunctionBeginUser;
28 3 : PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
29 3 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
30 3 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n1-D Laplacian, n=%" PetscInt_FMT "\n\n",n));
31 :
32 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
33 : Compute the operator matrix for the 1-D Laplacian
34 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
35 :
36 3 : PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
37 3 : PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n));
38 3 : PetscCall(MatSetFromOptions(A));
39 :
40 3 : PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
41 33 : for (i=Istart;i<Iend;i++) {
42 30 : if (i>0) PetscCall(MatSetValue(A,i,i-1,-1.0,INSERT_VALUES));
43 30 : if (i<n-1) PetscCall(MatSetValue(A,i,i+1,-1.0,INSERT_VALUES));
44 30 : PetscCall(MatSetValue(A,i,i,2.0,INSERT_VALUES));
45 : }
46 3 : PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
47 3 : PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
48 3 : PetscCall(MatCreateVecs(A,&v,&w));
49 3 : PetscCall(VecSet(v,1.0));
50 :
51 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
52 : Compute the split preconditioner matrix (one diagonal)
53 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
54 :
55 3 : PetscCall(MatCreate(PETSC_COMM_WORLD,&Pa));
56 3 : PetscCall(MatSetSizes(Pa,PETSC_DECIDE,PETSC_DECIDE,n,n));
57 3 : PetscCall(MatSetFromOptions(Pa));
58 :
59 3 : PetscCall(MatGetOwnershipRange(Pa,&Istart,&Iend));
60 33 : for (i=Istart;i<Iend;i++) PetscCall(MatSetValue(Pa,i,i,2.0,INSERT_VALUES));
61 3 : PetscCall(MatAssemblyBegin(Pa,MAT_FINAL_ASSEMBLY));
62 3 : PetscCall(MatAssemblyEnd(Pa,MAT_FINAL_ASSEMBLY));
63 :
64 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
65 : Create the spectral transformation object
66 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
67 :
68 3 : PetscCall(STCreate(PETSC_COMM_WORLD,&st));
69 3 : mat[0] = A;
70 3 : PetscCall(STSetMatrices(st,1,mat));
71 3 : mat[0] = Pa;
72 3 : PetscCall(STSetSplitPreconditioner(st,1,mat,SAME_NONZERO_PATTERN));
73 3 : PetscCall(STSetTransform(st,PETSC_TRUE));
74 3 : PetscCall(STSetFromOptions(st));
75 3 : PetscCall(STCayleySetAntishift(st,-0.3)); /* only relevant for cayley */
76 :
77 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
78 : Form the preconditioner matrix and print it
79 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
80 :
81 3 : PetscCall(PetscObjectTypeCompareAny((PetscObject)st,&flg,STSINVERT,STCAYLEY,""));
82 3 : if (flg) {
83 2 : PetscCall(STGetKSP(st,&ksp));
84 2 : PetscCall(KSPGetPC(ksp,&pc));
85 2 : PetscCall(STGetOperator(st,NULL));
86 2 : PetscCall(PCGetOperators(pc,NULL,&Pmat));
87 2 : PetscCall(MatView(Pmat,NULL));
88 : }
89 :
90 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
91 : Apply the operator
92 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
93 :
94 : /* sigma=0.0 */
95 3 : PetscCall(STSetUp(st));
96 3 : PetscCall(STGetType(st,&type));
97 3 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type));
98 3 : PetscCall(STApply(st,v,w));
99 3 : PetscCall(VecView(w,NULL));
100 :
101 : /* sigma=0.1 */
102 3 : sigma = 0.1;
103 3 : PetscCall(STSetShift(st,sigma));
104 3 : PetscCall(STGetShift(st,&sigma));
105 3 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma)));
106 3 : if (flg) {
107 2 : PetscCall(STGetOperator(st,NULL));
108 2 : PetscCall(PCGetOperators(pc,NULL,&Pmat));
109 2 : PetscCall(MatView(Pmat,NULL));
110 : }
111 3 : PetscCall(STApply(st,v,w));
112 3 : PetscCall(VecView(w,NULL));
113 :
114 3 : PetscCall(STDestroy(&st));
115 3 : PetscCall(MatDestroy(&A));
116 3 : PetscCall(MatDestroy(&Pa));
117 3 : PetscCall(VecDestroy(&v));
118 3 : PetscCall(VecDestroy(&w));
119 3 : PetscCall(SlepcFinalize());
120 : return 0;
121 : }
122 :
123 : /*TEST
124 :
125 : test:
126 : suffix: 1
127 : args: -st_type {{cayley shift sinvert}separate output}
128 : requires: !single
129 :
130 : TEST*/
|