Actual source code: test12.c
slepc-3.7.0 2016-05-16
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
8: SLEPc is free software: you can redistribute it and/or modify it under the
9: terms of version 3 of the GNU Lesser General Public License as published by
10: the Free Software Foundation.
12: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
13: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
14: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
15: more details.
17: You should have received a copy of the GNU Lesser General Public License
18: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
19: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
20: */
22: static char help[] = "Test block orthogonalization on a rank-deficient BV.\n\n";
24: #include <slepcbv.h>
28: int main(int argc,char **argv)
29: {
31: BV X,Z;
32: Mat M,R;
33: Vec v,w,t;
34: PetscInt i,j,n=20,k=8;
35: PetscViewer view;
36: PetscBool verbose;
37: PetscReal norm;
38: PetscScalar alpha;
40: SlepcInitialize(&argc,&argv,(char*)0,help);
41: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
42: PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL);
43: PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
44: PetscPrintf(PETSC_COMM_WORLD,"Test BV block orthogonalization (length %D, k=%D).\n",n,k);
45: if (k<6) SETERRQ(PETSC_COMM_SELF,1,"k must be at least 6");
47: /* Create template vector */
48: VecCreate(PETSC_COMM_WORLD,&t);
49: VecSetSizes(t,PETSC_DECIDE,n);
50: VecSetFromOptions(t);
52: /* Create BV object X */
53: BVCreate(PETSC_COMM_WORLD,&X);
54: PetscObjectSetName((PetscObject)X,"X");
55: BVSetSizesFromVec(X,t,k);
56: BVSetFromOptions(X);
58: /* Set up viewer */
59: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view);
60: if (verbose) {
61: PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB);
62: }
64: /* Fill X entries (first half) */
65: for (j=0;j<k/2;j++) {
66: BVGetColumn(X,j,&v);
67: VecSet(v,0.0);
68: for (i=0;i<=n/2;i++) {
69: if (i+j<n) {
70: alpha = (3.0*i+j-2)/(2*(i+j+1));
71: VecSetValue(v,i+j,alpha,INSERT_VALUES);
72: }
73: }
74: VecAssemblyBegin(v);
75: VecAssemblyEnd(v);
76: BVRestoreColumn(X,j,&v);
77: }
79: /* make middle column linearly dependent wrt columns 0 and 1 */
80: BVCopyColumn(X,0,j);
81: BVGetColumn(X,j,&v);
82: BVGetColumn(X,1,&w);
83: VecAXPY(v,0.5,w);
84: BVRestoreColumn(X,1,&w);
85: BVRestoreColumn(X,j,&v);
86: j++;
88: /* Fill X entries (second half) */
89: for (;j<k-1;j++) {
90: BVGetColumn(X,j,&v);
91: VecSet(v,0.0);
92: for (i=0;i<=n/2;i++) {
93: if (i+j<n) {
94: alpha = (3.0*i+j-2)/(2*(i+j+1));
95: VecSetValue(v,i+j,alpha,INSERT_VALUES);
96: }
97: }
98: VecAssemblyBegin(v);
99: VecAssemblyEnd(v);
100: BVRestoreColumn(X,j,&v);
101: }
103: /* make middle column linearly dependent wrt columns 1 and k/2+1 */
104: BVCopyColumn(X,1,j);
105: BVGetColumn(X,j,&v);
106: BVGetColumn(X,k/2+1,&w);
107: VecAXPY(v,-1.2,w);
108: BVRestoreColumn(X,k/2+1,&w);
109: BVRestoreColumn(X,j,&v);
111: if (verbose) {
112: BVView(X,view);
113: }
115: /* Create a copy on Z */
116: BVDuplicate(X,&Z);
117: PetscObjectSetName((PetscObject)Z,"Z");
118: BVCopy(X,Z);
120: /* Test BVOrthogonalize */
121: MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&R);
122: PetscObjectSetName((PetscObject)R,"R");
123: BVOrthogonalize(X,R);
124: if (verbose) {
125: BVView(X,view);
126: MatView(R,view);
127: }
129: /* Check orthogonality */
130: MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&M);
131: MatShift(M,1.0); /* set leading part to identity */
132: BVDot(X,X,M);
133: MatShift(M,-1.0);
134: MatNorm(M,NORM_1,&norm);
135: if (norm<100*PETSC_MACHINE_EPSILON) {
136: PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality < 100*eps\n");
137: } else {
138: PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g\n",(double)norm);
139: }
141: /* Check residual */
142: BVMult(Z,-1.0,1.0,X,R);
143: BVNorm(Z,NORM_FROBENIUS,&norm);
144: if (norm<100*PETSC_MACHINE_EPSILON) {
145: PetscPrintf(PETSC_COMM_WORLD,"Residual ||X-QR|| < 100*eps\n");
146: } else {
147: PetscPrintf(PETSC_COMM_WORLD,"Residual ||X-QR||: %g\n",(double)norm);
148: }
150: MatDestroy(&R);
151: MatDestroy(&M);
152: BVDestroy(&X);
153: BVDestroy(&Z);
154: VecDestroy(&t);
155: SlepcFinalize();
156: return ierr;
157: }