Actual source code: test12.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 block orthogonalization on a rank-deficient BV.\n\n";
13: #include <slepcbv.h>
15: int main(int argc,char **argv)
16: {
17: BV X,Z;
18: Mat M,R;
19: Vec v,w,t;
20: PetscInt i,j,n=20,k=8;
21: PetscViewer view;
22: PetscBool verbose;
23: PetscReal norm;
24: PetscScalar alpha;
26: SlepcInitialize(&argc,&argv,(char*)0,help);
27: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
28: PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL);
29: PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
30: PetscPrintf(PETSC_COMM_WORLD,"Test BV block orthogonalization (length %" PetscInt_FMT ", k=%" PetscInt_FMT ").\n",n,k);
33: /* Create template vector */
34: VecCreate(PETSC_COMM_WORLD,&t);
35: VecSetSizes(t,PETSC_DECIDE,n);
36: VecSetFromOptions(t);
38: /* Create BV object X */
39: BVCreate(PETSC_COMM_WORLD,&X);
40: PetscObjectSetName((PetscObject)X,"X");
41: BVSetSizesFromVec(X,t,k);
42: BVSetFromOptions(X);
44: /* Set up viewer */
45: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view);
46: if (verbose) PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB);
48: /* Fill X entries (first half) */
49: for (j=0;j<k/2;j++) {
50: BVGetColumn(X,j,&v);
51: VecSet(v,0.0);
52: for (i=0;i<=n/2;i++) {
53: if (i+j<n) {
54: alpha = (3.0*i+j-2)/(2*(i+j+1));
55: VecSetValue(v,i+j,alpha,INSERT_VALUES);
56: }
57: }
58: VecAssemblyBegin(v);
59: VecAssemblyEnd(v);
60: BVRestoreColumn(X,j,&v);
61: }
63: /* make middle column linearly dependent wrt columns 0 and 1 */
64: BVCopyColumn(X,0,j);
65: BVGetColumn(X,j,&v);
66: BVGetColumn(X,1,&w);
67: VecAXPY(v,0.5,w);
68: BVRestoreColumn(X,1,&w);
69: BVRestoreColumn(X,j,&v);
70: j++;
72: /* Fill X entries (second half) */
73: for (;j<k-1;j++) {
74: BVGetColumn(X,j,&v);
75: VecSet(v,0.0);
76: for (i=0;i<=n/2;i++) {
77: if (i+j<n) {
78: alpha = (3.0*i+j-2)/(2*(i+j+1));
79: VecSetValue(v,i+j,alpha,INSERT_VALUES);
80: }
81: }
82: VecAssemblyBegin(v);
83: VecAssemblyEnd(v);
84: BVRestoreColumn(X,j,&v);
85: }
87: /* make middle column linearly dependent wrt columns 1 and k/2+1 */
88: BVCopyColumn(X,1,j);
89: BVGetColumn(X,j,&v);
90: BVGetColumn(X,k/2+1,&w);
91: VecAXPY(v,-1.2,w);
92: BVRestoreColumn(X,k/2+1,&w);
93: BVRestoreColumn(X,j,&v);
95: if (verbose) BVView(X,view);
97: /* Create a copy on Z */
98: BVDuplicate(X,&Z);
99: PetscObjectSetName((PetscObject)Z,"Z");
100: BVCopy(X,Z);
102: /* Test BVOrthogonalize */
103: MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&R);
104: PetscObjectSetName((PetscObject)R,"R");
105: BVOrthogonalize(X,R);
106: if (verbose) {
107: BVView(X,view);
108: MatView(R,view);
109: }
111: /* Check orthogonality */
112: MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&M);
113: MatShift(M,1.0); /* set leading part to identity */
114: BVDot(X,X,M);
115: MatShift(M,-1.0);
116: MatNorm(M,NORM_1,&norm);
117: if (norm<100*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality < 100*eps\n");
118: else PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g\n",(double)norm);
120: /* Check residual */
121: BVMult(Z,-1.0,1.0,X,R);
122: BVNorm(Z,NORM_FROBENIUS,&norm);
123: if (norm<100*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"Residual ||X-QR|| < 100*eps\n");
124: else PetscPrintf(PETSC_COMM_WORLD,"Residual ||X-QR||: %g\n",(double)norm);
126: MatDestroy(&R);
127: MatDestroy(&M);
128: BVDestroy(&X);
129: BVDestroy(&Z);
130: VecDestroy(&t);
131: SlepcFinalize();
132: return 0;
133: }
135: /*TEST
137: test:
138: suffix: 1
139: nsize: 1
140: args: -bv_orthog_block gs -bv_type {{vecs contiguous svec mat}shared output}
141: output_file: output/test12_1.out
143: TEST*/