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