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 two matrices.\n\n";
12 :
13 : #include <slepcst.h>
14 :
15 3 : int main(int argc,char **argv)
16 : {
17 3 : Mat A,B,M,mat[2];
18 3 : ST st;
19 3 : Vec v,w;
20 3 : STType type;
21 3 : PetscScalar sigma,tau;
22 3 : PetscInt n=10,i,Istart,Iend;
23 :
24 3 : PetscFunctionBeginUser;
25 3 : PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
26 3 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
27 3 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n1-D Laplacian plus diagonal, n=%" PetscInt_FMT "\n\n",n));
28 :
29 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
30 : Compute the operator matrix for the 1-D Laplacian
31 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
32 :
33 3 : PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
34 3 : PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n));
35 3 : PetscCall(MatSetFromOptions(A));
36 :
37 3 : PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
38 3 : PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n));
39 3 : PetscCall(MatSetFromOptions(B));
40 :
41 3 : PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
42 33 : for (i=Istart;i<Iend;i++) {
43 30 : PetscCall(MatSetValue(A,i,i,2.0,INSERT_VALUES));
44 30 : if (i>0) {
45 27 : PetscCall(MatSetValue(A,i,i-1,-1.0,INSERT_VALUES));
46 27 : PetscCall(MatSetValue(B,i,i,(PetscScalar)i,INSERT_VALUES));
47 3 : } else PetscCall(MatSetValue(B,i,i,-1.0,INSERT_VALUES));
48 30 : if (i<n-1) PetscCall(MatSetValue(A,i,i+1,-1.0,INSERT_VALUES));
49 : }
50 3 : PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
51 3 : PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
52 3 : PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
53 3 : PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
54 3 : PetscCall(MatCreateVecs(A,&v,&w));
55 3 : PetscCall(VecSet(v,1.0));
56 :
57 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
58 : Create the spectral transformation object
59 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
60 :
61 3 : PetscCall(STCreate(PETSC_COMM_WORLD,&st));
62 3 : mat[0] = A;
63 3 : mat[1] = B;
64 3 : PetscCall(STSetMatrices(st,2,mat));
65 3 : PetscCall(STSetTransform(st,PETSC_TRUE));
66 3 : PetscCall(STSetFromOptions(st));
67 :
68 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
69 : Apply the transformed operator for several ST's
70 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
71 :
72 : /* shift, sigma=0.0 */
73 3 : PetscCall(STSetUp(st));
74 3 : PetscCall(STGetType(st,&type));
75 3 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type));
76 3 : PetscCall(STApply(st,v,w));
77 3 : PetscCall(VecView(w,NULL));
78 :
79 : /* shift, sigma=0.1 */
80 3 : sigma = 0.1;
81 3 : PetscCall(STSetShift(st,sigma));
82 3 : PetscCall(STGetShift(st,&sigma));
83 3 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma)));
84 3 : PetscCall(STApply(st,v,w));
85 3 : PetscCall(VecView(w,NULL));
86 :
87 : /* sinvert, sigma=0.1 */
88 3 : PetscCall(STPostSolve(st)); /* undo changes if inplace */
89 3 : PetscCall(STSetType(st,STSINVERT));
90 3 : PetscCall(STGetType(st,&type));
91 3 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type));
92 3 : PetscCall(STGetShift(st,&sigma));
93 3 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma)));
94 3 : PetscCall(STApply(st,v,w));
95 3 : PetscCall(VecView(w,NULL));
96 :
97 : /* sinvert, sigma=-0.5 */
98 3 : sigma = -0.5;
99 3 : PetscCall(STSetShift(st,sigma));
100 3 : PetscCall(STGetShift(st,&sigma));
101 3 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma)));
102 3 : PetscCall(STApply(st,v,w));
103 3 : PetscCall(VecView(w,NULL));
104 :
105 : /* cayley, sigma=-0.5, tau=-0.5 (equal to sigma by default) */
106 3 : PetscCall(STPostSolve(st)); /* undo changes if inplace */
107 3 : PetscCall(STSetType(st,STCAYLEY));
108 3 : PetscCall(STSetUp(st));
109 3 : PetscCall(STGetType(st,&type));
110 3 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type));
111 3 : PetscCall(STGetShift(st,&sigma));
112 3 : PetscCall(STCayleyGetAntishift(st,&tau));
113 3 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"With shift=%g, antishift=%g\n",(double)PetscRealPart(sigma),(double)PetscRealPart(tau)));
114 3 : PetscCall(STApply(st,v,w));
115 3 : PetscCall(VecView(w,NULL));
116 :
117 : /* cayley, sigma=1.1, tau=1.1 (still equal to sigma) */
118 3 : sigma = 1.1;
119 3 : PetscCall(STSetShift(st,sigma));
120 3 : PetscCall(STGetShift(st,&sigma));
121 3 : PetscCall(STCayleyGetAntishift(st,&tau));
122 3 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"With shift=%g, antishift=%g\n",(double)PetscRealPart(sigma),(double)PetscRealPart(tau)));
123 3 : PetscCall(STApply(st,v,w));
124 3 : PetscCall(VecView(w,NULL));
125 :
126 : /* cayley, sigma=1.1, tau=-1.0 */
127 3 : tau = -1.0;
128 3 : PetscCall(STCayleySetAntishift(st,tau));
129 3 : PetscCall(STGetShift(st,&sigma));
130 3 : PetscCall(STCayleyGetAntishift(st,&tau));
131 3 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"With shift=%g, antishift=%g\n",(double)PetscRealPart(sigma),(double)PetscRealPart(tau)));
132 3 : PetscCall(STApply(st,v,w));
133 3 : PetscCall(VecView(w,NULL));
134 :
135 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
136 : Check inner product matrix in Cayley
137 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
138 3 : PetscCall(STGetBilinearForm(st,&M));
139 3 : PetscCall(MatMult(M,v,w));
140 3 : PetscCall(VecView(w,NULL));
141 :
142 3 : PetscCall(STDestroy(&st));
143 3 : PetscCall(MatDestroy(&A));
144 3 : PetscCall(MatDestroy(&B));
145 3 : PetscCall(MatDestroy(&M));
146 3 : PetscCall(VecDestroy(&v));
147 3 : PetscCall(VecDestroy(&w));
148 3 : PetscCall(SlepcFinalize());
149 : return 0;
150 : }
151 :
152 : /*TEST
153 :
154 : test:
155 : suffix: 1
156 : args: -st_matmode {{copy inplace shell}}
157 : requires: !single
158 :
159 : TEST*/
|