Actual source code: test11.c
slepc-3.19.0 2023-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;
24: PetscFunctionBeginUser;
25: PetscCall(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: PetscCall(MatDenseRestoreArrayRead(A,&pA));
35: PetscFunctionReturn(PETSC_SUCCESS);
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: PetscFunctionBeginUser;
50: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
51: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
52: PetscCall(PetscOptionsGetInt(NULL,NULL,"-l",&l,NULL));
53: PetscCall(PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL));
54: PetscCall(PetscOptionsHasName(NULL,NULL,"-withb",&withb));
55: PetscCall(PetscOptionsHasName(NULL,NULL,"-resid",&resid));
56: PetscCall(PetscOptionsHasName(NULL,NULL,"-rand",&rand));
57: PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose));
58: PetscCall(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":""));
60: /* Create template vector */
61: PetscCall(VecCreate(PETSC_COMM_WORLD,&t));
62: PetscCall(VecSetSizes(t,PETSC_DECIDE,n));
63: PetscCall(VecSetFromOptions(t));
65: /* Create BV object X */
66: PetscCall(BVCreate(PETSC_COMM_WORLD,&X));
67: PetscCall(PetscObjectSetName((PetscObject)X,"X"));
68: PetscCall(BVSetSizesFromVec(X,t,k));
69: PetscCall(BVSetFromOptions(X));
71: /* Set up viewer */
72: PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view));
73: if (verbose) PetscCall(PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB));
75: /* Fill X entries */
76: if (rand) PetscCall(BVSetRandom(X));
77: else {
78: for (j=0;j<k;j++) {
79: PetscCall(BVGetColumn(X,j,&v));
80: PetscCall(VecSet(v,0.0));
81: for (i=0;i<=n/2;i++) {
82: if (i+j<n) {
83: alpha = (3.0*i+j-2)/(2*(i+j+1));
84: PetscCall(VecSetValue(v,i+j,alpha,INSERT_VALUES));
85: }
86: }
87: PetscCall(VecAssemblyBegin(v));
88: PetscCall(VecAssemblyEnd(v));
89: PetscCall(BVRestoreColumn(X,j,&v));
90: }
91: }
92: if (verbose) PetscCall(BVView(X,view));
94: if (withb) {
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));
111: PetscCall(BVSetMatrix(X,B,PETSC_FALSE));
112: }
114: /* Create copy on Y */
115: PetscCall(BVDuplicate(X,&Y));
116: PetscCall(PetscObjectSetName((PetscObject)Y,"Y"));
117: PetscCall(BVCopy(X,Y));
118: PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&M));
120: if (resid) {
121: /* Create matrix R to store triangular factor */
122: PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&R));
123: PetscCall(PetscObjectSetName((PetscObject)R,"R"));
124: }
126: if (l>0) {
127: /* First orthogonalize leading columns */
128: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Orthogonalizing leading columns\n"));
129: PetscCall(BVSetActiveColumns(Y,0,l));
130: PetscCall(BVSetActiveColumns(X,0,l));
131: PetscCall(BVOrthogonalize(Y,R));
132: if (verbose) {
133: PetscCall(BVView(Y,view));
134: if (resid) PetscCall(MatView(R,view));
135: }
137: if (withb) {
138: /* Extract cached BV and check it is equal to B*X */
139: PetscCall(BVGetCachedBV(Y,&cached));
140: PetscCall(BVDuplicate(X,&Z));
141: PetscCall(BVSetMatrix(Z,NULL,PETSC_FALSE));
142: PetscCall(BVSetActiveColumns(Z,0,l));
143: PetscCall(BVCopy(X,Z));
144: PetscCall(BVMatMult(X,B,Z));
145: PetscCall(BVMult(Z,-1.0,1.0,cached,NULL));
146: PetscCall(BVNorm(Z,NORM_FROBENIUS,&norm));
147: if (norm<100*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD," Difference ||cached-BX|| < 100*eps\n"));
148: else PetscCall(PetscPrintf(PETSC_COMM_WORLD," Difference ||cached-BX||: %g\n",(double)norm));
149: PetscCall(BVDestroy(&Z));
150: }
152: /* Check orthogonality */
153: PetscCall(BVDot(Y,Y,M));
154: PetscCall(MyMatNorm(M,k,0,l,1.0,&norm));
155: if (norm<100*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD," Level of orthogonality of Q1 < 100*eps\n"));
156: else PetscCall(PetscPrintf(PETSC_COMM_WORLD," Level of orthogonality of Q1: %g\n",(double)norm));
158: if (resid) {
159: /* Check residual */
160: PetscCall(BVDuplicate(X,&Z));
161: PetscCall(BVSetMatrix(Z,NULL,PETSC_FALSE));
162: PetscCall(BVSetActiveColumns(Z,0,l));
163: PetscCall(BVCopy(X,Z));
164: PetscCall(BVMult(Z,-1.0,1.0,Y,R));
165: PetscCall(BVNorm(Z,NORM_FROBENIUS,&norm));
166: if (norm<100*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD," Residual ||X1-Q1*R11|| < 100*eps\n"));
167: else PetscCall(PetscPrintf(PETSC_COMM_WORLD," Residual ||X1-Q1*R11||: %g\n",(double)norm));
168: PetscCall(BVDestroy(&Z));
169: }
171: }
173: /* Now orthogonalize the rest of columns */
174: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Orthogonalizing active columns\n"));
175: PetscCall(BVSetActiveColumns(Y,l,k));
176: PetscCall(BVSetActiveColumns(X,l,k));
177: PetscCall(BVOrthogonalize(Y,R));
178: if (verbose) {
179: PetscCall(BVView(Y,view));
180: if (resid) PetscCall(MatView(R,view));
181: }
183: if (l>0) {
184: /* Check orthogonality */
185: PetscCall(BVDot(Y,Y,M));
186: PetscCall(MyMatNorm(M,k,l,k,1.0,&norm));
187: if (norm<100*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD," Level of orthogonality of Q2 < 100*eps\n"));
188: else PetscCall(PetscPrintf(PETSC_COMM_WORLD," Level of orthogonality of Q2: %g\n",(double)norm));
189: }
191: /* Check the complete decomposition */
192: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Overall decomposition\n"));
193: PetscCall(BVSetActiveColumns(Y,0,k));
194: PetscCall(BVSetActiveColumns(X,0,k));
196: /* Check orthogonality */
197: PetscCall(BVDot(Y,Y,M));
198: PetscCall(MyMatNorm(M,k,0,k,1.0,&norm));
199: if (norm<100*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD," Level of orthogonality of Q < 100*eps\n"));
200: else PetscCall(PetscPrintf(PETSC_COMM_WORLD," Level of orthogonality of Q: %g\n",(double)norm));
202: if (resid) {
203: /* Check residual */
204: PetscCall(BVMult(X,-1.0,1.0,Y,R));
205: PetscCall(BVSetMatrix(X,NULL,PETSC_FALSE));
206: PetscCall(BVNorm(X,NORM_FROBENIUS,&norm));
207: if (norm<100*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD," Residual ||X-Q*R|| < 100*eps\n"));
208: else PetscCall(PetscPrintf(PETSC_COMM_WORLD," Residual ||X-Q*R||: %g\n",(double)norm));
209: PetscCall(MatDestroy(&R));
210: }
212: if (B) PetscCall(MatDestroy(&B));
213: PetscCall(MatDestroy(&M));
214: PetscCall(BVDestroy(&X));
215: PetscCall(BVDestroy(&Y));
216: PetscCall(VecDestroy(&t));
217: PetscCall(SlepcFinalize());
218: return 0;
219: }
221: /*TEST
223: testset:
224: args: -bv_orthog_block {{gs chol tsqr tsqrchol svqb}}
225: nsize: 2
226: output_file: output/test11_1.out
227: test:
228: suffix: 1
229: args: -bv_type {{vecs contiguous svec mat}shared output}
230: test:
231: suffix: 1_cuda
232: args: -bv_type svec -vec_type cuda
233: requires: cuda
235: testset:
236: args: -withb -bv_orthog_block {{gs chol svqb}}
237: nsize: 2
238: output_file: output/test11_4.out
239: test:
240: suffix: 4
241: args: -bv_type {{vecs contiguous svec mat}shared output}
242: test:
243: suffix: 4_cuda
244: args: -bv_type svec -vec_type cuda -mat_type aijcusparse
245: requires: cuda
247: testset:
248: args: -resid -bv_orthog_block {{gs chol tsqr tsqrchol svqb}}
249: nsize: 2
250: output_file: output/test11_6.out
251: test:
252: suffix: 6
253: args: -bv_type {{vecs contiguous svec mat}shared output}
254: test:
255: suffix: 6_cuda
256: args: -bv_type svec -vec_type cuda
257: requires: cuda
259: testset:
260: args: -resid -withb -bv_orthog_block {{gs chol svqb}}
261: nsize: 2
262: output_file: output/test11_9.out
263: test:
264: suffix: 9
265: args: -bv_type {{vecs contiguous svec mat}shared output}
266: test:
267: suffix: 9_cuda
268: args: -bv_type svec -vec_type cuda -mat_type aijcusparse
269: requires: cuda
271: testset:
272: args: -bv_orthog_block tsqr
273: nsize: 7
274: output_file: output/test11_1.out
275: test:
276: suffix: 11
277: args: -bv_type {{vecs contiguous svec mat}shared output}
278: requires: !defined(PETSCTEST_VALGRIND)
279: test:
280: suffix: 11_cuda
281: TODO: too many processes accessing the GPU
282: args: -bv_type svec -vec_type cuda
283: requires: cuda !defined(PETSCTEST_VALGRIND)
285: testset:
286: args: -resid -n 180 -l 0 -k 7 -bv_orthog_block tsqr
287: nsize: 9
288: output_file: output/test11_12.out
289: test:
290: suffix: 12
291: args: -bv_type {{vecs contiguous svec mat}shared output}
292: requires: !single !defined(PETSCTEST_VALGRIND)
293: test:
294: suffix: 12_cuda
295: TODO: too many processes accessing the GPU
296: args: -bv_type svec -vec_type cuda
297: requires: cuda !single !defined(PETSCTEST_VALGRIND)
299: TEST*/