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 BVNormalize().\n\n";
12 :
13 : #include <slepcbv.h>
14 :
15 12 : int main(int argc,char **argv)
16 : {
17 12 : BV X,Y,Z;
18 12 : Mat B;
19 12 : Vec v,t;
20 12 : PetscInt i,j,n=20,k=8,l=3,Istart,Iend;
21 12 : PetscViewer view;
22 12 : PetscBool verbose;
23 12 : PetscReal norm,error;
24 12 : PetscScalar alpha;
25 : #if !defined(PETSC_USE_COMPLEX)
26 12 : PetscScalar *eigi;
27 12 : PetscRandom rand;
28 12 : PetscReal normr,normi;
29 12 : Vec vi;
30 : #endif
31 :
32 12 : PetscFunctionBeginUser;
33 12 : PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
34 12 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
35 12 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL));
36 12 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-l",&l,NULL));
37 12 : PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose));
38 12 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test BV normalization with %" PetscInt_FMT " columns of length %" PetscInt_FMT ".\n",k,n));
39 :
40 : /* Create template vector */
41 12 : PetscCall(VecCreate(PETSC_COMM_WORLD,&t));
42 12 : PetscCall(VecSetSizes(t,PETSC_DECIDE,n));
43 12 : PetscCall(VecSetFromOptions(t));
44 :
45 : /* Create BV object X */
46 12 : PetscCall(BVCreate(PETSC_COMM_WORLD,&X));
47 12 : PetscCall(PetscObjectSetName((PetscObject)X,"X"));
48 12 : PetscCall(BVSetSizesFromVec(X,t,k));
49 12 : PetscCall(BVSetFromOptions(X));
50 :
51 : /* Set up viewer */
52 12 : PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view));
53 12 : if (verbose) PetscCall(PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB));
54 :
55 : /* Fill X entries */
56 192 : for (j=0;j<k;j++) {
57 180 : PetscCall(BVGetColumn(X,j,&v));
58 180 : PetscCall(VecSet(v,0.0));
59 22860 : for (i=0;i<=n/2;i++) {
60 22680 : if (i+j<n) {
61 22680 : alpha = (3.0*i+j-2)/(2*(i+j+1));
62 22680 : PetscCall(VecSetValue(v,i+j,alpha,INSERT_VALUES));
63 : }
64 : }
65 180 : PetscCall(VecAssemblyBegin(v));
66 180 : PetscCall(VecAssemblyEnd(v));
67 180 : PetscCall(BVRestoreColumn(X,j,&v));
68 : }
69 12 : if (verbose) PetscCall(BVView(X,view));
70 :
71 : /* Create copies on Y and Z */
72 12 : PetscCall(BVDuplicate(X,&Y));
73 12 : PetscCall(PetscObjectSetName((PetscObject)Y,"Y"));
74 12 : PetscCall(BVCopy(X,Y));
75 12 : PetscCall(BVDuplicate(X,&Z));
76 12 : PetscCall(PetscObjectSetName((PetscObject)Z,"Z"));
77 12 : PetscCall(BVCopy(X,Z));
78 12 : PetscCall(BVSetActiveColumns(X,l,k));
79 12 : PetscCall(BVSetActiveColumns(Y,l,k));
80 12 : PetscCall(BVSetActiveColumns(Z,l,k));
81 :
82 : /* Test BVNormalize */
83 12 : PetscCall(BVNormalize(X,NULL));
84 12 : if (verbose) PetscCall(BVView(X,view));
85 :
86 : /* Check unit norm of columns */
87 12 : error = 0.0;
88 120 : for (j=l;j<k;j++) {
89 108 : PetscCall(BVNormColumn(X,j,NORM_2,&norm));
90 129 : error = PetscMax(error,PetscAbsReal(norm-PetscRealConstant(1.0)));
91 : }
92 12 : if (error<100*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Deviation from normalized vectors < 100*eps\n"));
93 0 : else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Deviation from normalized vectors: %g\n",(double)norm));
94 :
95 : /* Create inner product matrix */
96 12 : PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
97 12 : PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n));
98 12 : PetscCall(MatSetFromOptions(B));
99 12 : PetscCall(PetscObjectSetName((PetscObject)B,"B"));
100 :
101 12 : PetscCall(MatGetOwnershipRange(B,&Istart,&Iend));
102 2012 : for (i=Istart;i<Iend;i++) {
103 2000 : if (i>0) PetscCall(MatSetValue(B,i,i-1,-1.0,INSERT_VALUES));
104 2000 : if (i<n-1) PetscCall(MatSetValue(B,i,i+1,-1.0,INSERT_VALUES));
105 2000 : PetscCall(MatSetValue(B,i,i,2.0,INSERT_VALUES));
106 : }
107 12 : PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
108 12 : PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
109 12 : if (verbose) PetscCall(MatView(B,view));
110 :
111 : /* Test BVNormalize with B-norm */
112 12 : PetscCall(BVSetMatrix(Y,B,PETSC_FALSE));
113 12 : PetscCall(BVNormalize(Y,NULL));
114 12 : if (verbose) PetscCall(BVView(Y,view));
115 :
116 : /* Check unit B-norm of columns */
117 12 : error = 0.0;
118 120 : for (j=l;j<k;j++) {
119 108 : PetscCall(BVNormColumn(Y,j,NORM_2,&norm));
120 128 : error = PetscMax(error,PetscAbsReal(norm-PetscRealConstant(1.0)));
121 : }
122 12 : if (error<100*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Deviation from B-normalized vectors < 100*eps\n"));
123 0 : else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Deviation from B-normalized vectors: %g\n",(double)norm));
124 :
125 : #if !defined(PETSC_USE_COMPLEX)
126 : /* fill imaginary parts */
127 12 : PetscCall(PetscCalloc1(k,&eigi));
128 12 : PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&rand));
129 12 : PetscCall(PetscRandomSetFromOptions(rand));
130 36 : for (j=l+1;j<k-1;j+=5) {
131 24 : PetscCall(PetscRandomGetValue(rand,&alpha));
132 24 : eigi[j] = alpha;
133 24 : eigi[j+1] = -alpha;
134 : }
135 12 : PetscCall(PetscRandomDestroy(&rand));
136 12 : if (verbose) {
137 0 : PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,1,k,eigi,&v));
138 0 : PetscCall(VecView(v,view));
139 0 : PetscCall(VecDestroy(&v));
140 : }
141 :
142 : /* Test BVNormalize with complex conjugate columns */
143 12 : PetscCall(BVNormalize(Z,eigi));
144 12 : if (verbose) PetscCall(BVView(Z,view));
145 :
146 : /* Check unit norm of (complex conjugate) columns */
147 12 : error = 0.0;
148 96 : for (j=l;j<k;j++) {
149 84 : if (eigi[j]) {
150 24 : PetscCall(BVGetColumn(Z,j,&v));
151 24 : PetscCall(BVGetColumn(Z,j+1,&vi));
152 24 : PetscCall(VecNormBegin(v,NORM_2,&normr));
153 24 : PetscCall(VecNormBegin(vi,NORM_2,&normi));
154 24 : PetscCall(VecNormEnd(v,NORM_2,&normr));
155 24 : PetscCall(VecNormEnd(vi,NORM_2,&normi));
156 24 : PetscCall(BVRestoreColumn(Z,j+1,&vi));
157 24 : PetscCall(BVRestoreColumn(Z,j,&v));
158 24 : norm = SlepcAbsEigenvalue(normr,normi);
159 24 : j++;
160 : } else {
161 60 : PetscCall(BVGetColumn(Z,j,&v));
162 60 : PetscCall(VecNorm(v,NORM_2,&norm));
163 60 : PetscCall(BVRestoreColumn(Z,j,&v));
164 : }
165 96 : error = PetscMax(error,PetscAbsReal(norm-PetscRealConstant(1.0)));
166 : }
167 12 : if (error<100*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Deviation from normalized conjugate vectors < 100*eps\n"));
168 0 : else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Deviation from normalized conjugate vectors: %g\n",(double)norm));
169 12 : PetscCall(PetscFree(eigi));
170 : #endif
171 :
172 12 : PetscCall(BVDestroy(&X));
173 12 : PetscCall(BVDestroy(&Y));
174 12 : PetscCall(BVDestroy(&Z));
175 12 : PetscCall(MatDestroy(&B));
176 12 : PetscCall(VecDestroy(&t));
177 12 : PetscCall(SlepcFinalize());
178 : return 0;
179 : }
180 :
181 : /*TEST
182 :
183 : testset:
184 : args: -n 250 -l 6 -k 15
185 : nsize: {{1 2}}
186 : requires: !complex
187 : output_file: output/test18_1.out
188 : test:
189 : suffix: 1
190 : args: -bv_type {{vecs contiguous svec mat}}
191 : test:
192 : suffix: 1_cuda
193 : args: -bv_type {{svec mat}} -vec_type cuda
194 : requires: cuda
195 : test:
196 : suffix: 1_hip
197 : args: -bv_type {{svec mat}} -vec_type hip
198 : requires: hip
199 :
200 : testset:
201 : args: -n 250 -l 6 -k 15
202 : nsize: {{1 2}}
203 : requires: complex
204 : output_file: output/test18_1_complex.out
205 : test:
206 : suffix: 1_complex
207 : args: -bv_type {{vecs contiguous svec mat}}
208 : test:
209 : suffix: 1_cuda_complex
210 : args: -bv_type {{svec mat}} -vec_type cuda
211 : requires: cuda
212 : test:
213 : suffix: 1_hip_complex
214 : args: -bv_type {{svec mat}} -vec_type hip
215 : requires: hip
216 :
217 : TEST*/
|