Actual source code: test11.c
slepc-3.17.0 2022-03-31
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 BV block orthogonalization.\n\n";
13: #include <slepcbv.h>
15: /*
16: Compute the Frobenius norm ||A(l:k,l:k)-diag||_F
17: */
18: PetscErrorCode MyMatNorm(Mat A,PetscInt lda,PetscInt l,PetscInt k,PetscScalar diag,PetscReal *norm)
19: {
20: PetscInt i,j;
21: const PetscScalar *pA;
22: PetscReal s,val;
25: MatDenseGetArrayRead(A,&pA);
26: s = 0.0;
27: for (i=l;i<k;i++) {
28: for (j=l;j<k;j++) {
29: val = (i==j)? PetscAbsScalar(pA[i+j*lda]-diag): PetscAbsScalar(pA[i+j*lda]);
30: s += val*val;
31: }
32: }
33: *norm = PetscSqrtReal(s);
34: MatDenseRestoreArrayRead(A,&pA);
35: PetscFunctionReturn(0);
36: }
38: int main(int argc,char **argv)
39: {
40: BV X,Y,Z,cached;
41: Mat B=NULL,M,R=NULL;
42: Vec v,t;
43: PetscInt i,j,n=20,l=2,k=8,Istart,Iend;
44: PetscViewer view;
45: PetscBool withb,resid,rand,verbose;
46: PetscReal norm;
47: PetscScalar alpha;
49: SlepcInitialize(&argc,&argv,(char*)0,help);
50: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
51: PetscOptionsGetInt(NULL,NULL,"-l",&l,NULL);
52: PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL);
53: PetscOptionsHasName(NULL,NULL,"-withb",&withb);
54: PetscOptionsHasName(NULL,NULL,"-resid",&resid);
55: PetscOptionsHasName(NULL,NULL,"-rand",&rand);
56: PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
57: PetscPrintf(PETSC_COMM_WORLD,"Test BV block orthogonalization (length %" PetscInt_FMT ", l=%" PetscInt_FMT ", k=%" PetscInt_FMT ")%s.\n",n,l,k,withb?" with non-standard inner product":"");
59: /* Create template vector */
60: VecCreate(PETSC_COMM_WORLD,&t);
61: VecSetSizes(t,PETSC_DECIDE,n);
62: VecSetFromOptions(t);
64: /* Create BV object X */
65: BVCreate(PETSC_COMM_WORLD,&X);
66: PetscObjectSetName((PetscObject)X,"X");
67: BVSetSizesFromVec(X,t,k);
68: BVSetFromOptions(X);
70: /* Set up viewer */
71: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view);
72: if (verbose) PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB);
74: /* Fill X entries */
75: if (rand) BVSetRandom(X);
76: else {
77: for (j=0;j<k;j++) {
78: BVGetColumn(X,j,&v);
79: VecSet(v,0.0);
80: for (i=0;i<=n/2;i++) {
81: if (i+j<n) {
82: alpha = (3.0*i+j-2)/(2*(i+j+1));
83: VecSetValue(v,i+j,alpha,INSERT_VALUES);
84: }
85: }
86: VecAssemblyBegin(v);
87: VecAssemblyEnd(v);
88: BVRestoreColumn(X,j,&v);
89: }
90: }
91: if (verbose) BVView(X,view);
93: if (withb) {
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);
110: BVSetMatrix(X,B,PETSC_FALSE);
111: }
113: /* Create copy on Y */
114: BVDuplicate(X,&Y);
115: PetscObjectSetName((PetscObject)Y,"Y");
116: BVCopy(X,Y);
117: MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&M);
119: if (resid) {
120: /* Create matrix R to store triangular factor */
121: MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&R);
122: PetscObjectSetName((PetscObject)R,"R");
123: }
125: if (l>0) {
126: /* First orthogonalize leading columns */
127: PetscPrintf(PETSC_COMM_WORLD,"Orthogonalizing leading columns\n");
128: BVSetActiveColumns(Y,0,l);
129: BVSetActiveColumns(X,0,l);
130: BVOrthogonalize(Y,R);
131: if (verbose) {
132: BVView(Y,view);
133: if (resid) MatView(R,view);
134: }
136: if (withb) {
137: /* Extract cached BV and check it is equal to B*X */
138: BVGetCachedBV(Y,&cached);
139: BVDuplicate(X,&Z);
140: BVSetMatrix(Z,NULL,PETSC_FALSE);
141: BVSetActiveColumns(Z,0,l);
142: BVCopy(X,Z);
143: BVMatMult(X,B,Z);
144: BVMult(Z,-1.0,1.0,cached,NULL);
145: BVNorm(Z,NORM_FROBENIUS,&norm);
146: if (norm<100*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD," Difference ||cached-BX|| < 100*eps\n");
147: else PetscPrintf(PETSC_COMM_WORLD," Difference ||cached-BX||: %g\n",(double)norm);
148: BVDestroy(&Z);
149: }
151: /* Check orthogonality */
152: BVDot(Y,Y,M);
153: MyMatNorm(M,k,0,l,1.0,&norm);
154: if (norm<100*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD," Level of orthogonality of Q1 < 100*eps\n");
155: else PetscPrintf(PETSC_COMM_WORLD," Level of orthogonality of Q1: %g\n",(double)norm);
157: if (resid) {
158: /* Check residual */
159: BVDuplicate(X,&Z);
160: BVSetMatrix(Z,NULL,PETSC_FALSE);
161: BVSetActiveColumns(Z,0,l);
162: BVCopy(X,Z);
163: BVMult(Z,-1.0,1.0,Y,R);
164: BVNorm(Z,NORM_FROBENIUS,&norm);
165: if (norm<100*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD," Residual ||X1-Q1*R11|| < 100*eps\n");
166: else PetscPrintf(PETSC_COMM_WORLD," Residual ||X1-Q1*R11||: %g\n",(double)norm);
167: BVDestroy(&Z);
168: }
170: }
172: /* Now orthogonalize the rest of columns */
173: PetscPrintf(PETSC_COMM_WORLD,"Orthogonalizing active columns\n");
174: BVSetActiveColumns(Y,l,k);
175: BVSetActiveColumns(X,l,k);
176: BVOrthogonalize(Y,R);
177: if (verbose) {
178: BVView(Y,view);
179: if (resid) MatView(R,view);
180: }
182: if (l>0) {
183: /* Check orthogonality */
184: BVDot(Y,Y,M);
185: MyMatNorm(M,k,l,k,1.0,&norm);
186: if (norm<100*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD," Level of orthogonality of Q2 < 100*eps\n");
187: else PetscPrintf(PETSC_COMM_WORLD," Level of orthogonality of Q2: %g\n",(double)norm);
188: }
190: /* Check the complete decomposition */
191: PetscPrintf(PETSC_COMM_WORLD,"Overall decomposition\n");
192: BVSetActiveColumns(Y,0,k);
193: BVSetActiveColumns(X,0,k);
195: /* Check orthogonality */
196: BVDot(Y,Y,M);
197: MyMatNorm(M,k,0,k,1.0,&norm);
198: if (norm<100*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD," Level of orthogonality of Q < 100*eps\n");
199: else PetscPrintf(PETSC_COMM_WORLD," Level of orthogonality of Q: %g\n",(double)norm);
201: if (resid) {
202: /* Check residual */
203: BVMult(X,-1.0,1.0,Y,R);
204: BVSetMatrix(X,NULL,PETSC_FALSE);
205: BVNorm(X,NORM_FROBENIUS,&norm);
206: if (norm<100*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD," Residual ||X-Q*R|| < 100*eps\n");
207: else PetscPrintf(PETSC_COMM_WORLD," Residual ||X-Q*R||: %g\n",(double)norm);
208: MatDestroy(&R);
209: }
211: if (B) MatDestroy(&B);
212: MatDestroy(&M);
213: BVDestroy(&X);
214: BVDestroy(&Y);
215: VecDestroy(&t);
216: SlepcFinalize();
217: return 0;
218: }
220: /*TEST
222: testset:
223: args: -bv_orthog_block {{gs chol tsqr tsqrchol svqb}}
224: nsize: 2
225: output_file: output/test11_1.out
226: test:
227: suffix: 1
228: args: -bv_type {{vecs contiguous svec mat}shared output}
229: test:
230: suffix: 1_cuda
231: args: -bv_type svec -vec_type cuda
232: requires: cuda
234: testset:
235: args: -withb -bv_orthog_block {{gs chol svqb}}
236: nsize: 2
237: output_file: output/test11_4.out
238: test:
239: suffix: 4
240: args: -bv_type {{vecs contiguous svec mat}shared output}
241: test:
242: suffix: 4_cuda
243: args: -bv_type svec -vec_type cuda -mat_type aijcusparse
244: requires: cuda
246: testset:
247: args: -resid -bv_orthog_block {{gs chol tsqr tsqrchol svqb}}
248: nsize: 2
249: output_file: output/test11_6.out
250: test:
251: suffix: 6
252: args: -bv_type {{vecs contiguous svec mat}shared output}
253: test:
254: suffix: 6_cuda
255: args: -bv_type svec -vec_type cuda
256: requires: cuda
258: testset:
259: args: -resid -withb -bv_orthog_block {{gs chol svqb}}
260: nsize: 2
261: output_file: output/test11_9.out
262: test:
263: suffix: 9
264: args: -bv_type {{vecs contiguous svec mat}shared output}
265: test:
266: suffix: 9_cuda
267: args: -bv_type svec -vec_type cuda -mat_type aijcusparse
268: requires: cuda
270: testset:
271: args: -bv_orthog_block tsqr
272: nsize: 7
273: output_file: output/test11_1.out
274: test:
275: suffix: 11
276: args: -bv_type {{vecs contiguous svec mat}shared output}
277: requires: !valgrind
278: test:
279: suffix: 11_cuda
280: args: -bv_type svec -vec_type cuda
281: requires: cuda !valgrind
283: testset:
284: args: -resid -n 180 -l 0 -k 7 -bv_orthog_block tsqr
285: nsize: 9
286: output_file: output/test11_12.out
287: test:
288: suffix: 12
289: args: -bv_type {{vecs contiguous svec mat}shared output}
290: requires: !single !valgrind
291: test:
292: suffix: 12_cuda
293: args: -bv_type svec -vec_type cuda
294: requires: !single !valgrind cuda
296: TEST*/