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