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