Actual source code: test5.c
slepc-3.6.3 2016-03-29
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2015, 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 BV operations with indefinite inner product.\n\n";
24: #include <slepcbv.h>
28: int main(int argc,char **argv)
29: {
31: Vec t,v;
32: Mat B,M;
33: Vec omega;
34: BV X,Y;
35: PetscInt i,j,n=10,k=5,Istart,Iend;
36: PetscScalar alpha;
37: PetscReal nrm;
38: PetscViewer view;
39: PetscBool verbose;
41: SlepcInitialize(&argc,&argv,(char*)0,help);
42: PetscOptionsGetInt(NULL,"-n",&n,NULL);
43: PetscOptionsGetInt(NULL,"-k",&k,NULL);
44: PetscOptionsHasName(NULL,"-verbose",&verbose);
45: PetscPrintf(PETSC_COMM_WORLD,"Test BV with indefinite inner product (n=%D, k=%D).\n",n,k);
47: /* Create inner product matrix (standard involutionary permutation) */
48: MatCreate(PETSC_COMM_WORLD,&B);
49: MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n);
50: MatSetFromOptions(B);
51: MatSetUp(B);
52: PetscObjectSetName((PetscObject)B,"B");
54: MatGetOwnershipRange(B,&Istart,&Iend);
55: for (i=Istart;i<Iend;i++) {
56: MatSetValue(B,i,n-i-1,1.0,INSERT_VALUES);
57: }
58: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
59: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
60: MatCreateVecs(B,&t,NULL);
62: /* Create BV object X */
63: BVCreate(PETSC_COMM_WORLD,&X);
64: PetscObjectSetName((PetscObject)X,"X");
65: BVSetSizesFromVec(X,t,k);
66: BVSetFromOptions(X);
67: BVSetMatrix(X,B,PETSC_TRUE);
69: /* Set up viewer */
70: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view);
71: if (verbose) {
72: PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB);
73: }
75: /* Fill X entries */
76: for (j=0;j<k;j++) {
77: BVGetColumn(X,j,&v);
78: VecSet(v,-1.0);
79: for (i=0;i<4;i++) {
80: if (i+j<n) {
81: VecSetValue(v,i+j,(PetscScalar)(3*i+j-2),INSERT_VALUES);
82: }
83: }
84: VecAssemblyBegin(v);
85: VecAssemblyEnd(v);
86: BVRestoreColumn(X,j,&v);
87: }
88: if (verbose) {
89: MatView(B,view);
90: BVView(X,view);
91: }
93: /* Test BVNormColumn */
94: BVNormColumn(X,0,NORM_2,&nrm);
95: PetscPrintf(PETSC_COMM_WORLD,"B-Norm or X[0] = %g\n",(double)nrm);
97: /* Test BVOrthogonalizeColumn */
98: for (j=0;j<k;j++) {
99: BVOrthogonalizeColumn(X,j,NULL,&nrm,NULL);
100: alpha = 1.0/nrm;
101: BVScaleColumn(X,j,alpha);
102: }
103: if (verbose) {
104: BVView(X,view);
105: }
107: /* Create a copy on Y */
108: BVDuplicate(X,&Y);
109: PetscObjectSetName((PetscObject)Y,"Y");
110: BVCopy(X,Y);
112: /* Check orthogonality */
113: MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&M);
114: BVDot(Y,Y,M);
115: VecCreateSeq(PETSC_COMM_SELF,k,&omega);
116: BVGetSignature(Y,omega);
117: VecScale(omega,-1.0);
118: MatDiagonalSet(M,omega,ADD_VALUES);
119: VecDestroy(&omega);
120: MatNorm(M,NORM_1,&nrm);
121: if (nrm<100*PETSC_MACHINE_EPSILON) {
122: PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality < 100*eps\n");
123: } else {
124: PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g\n",(double)nrm);
125: }
127: BVDestroy(&X);
128: BVDestroy(&Y);
129: MatDestroy(&M);
130: MatDestroy(&B);
131: VecDestroy(&t);
132: SlepcFinalize();
133: return 0;
134: }