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