Actual source code: test9.c
slepc-3.5.2 2014-10-10
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2014, 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 matrix projection.\n\n";
24: #include <slepcbv.h>
28: int main(int argc,char **argv)
29: {
31: Vec t,v;
32: Mat B,G,H0,H1;
33: BV X,Y,Z;
34: PetscInt i,j,n=20,kx=6,lx=3,ky=5,ly=2,Istart,Iend,col[5];
35: PetscScalar value[] = { -1, 1, 1, 1, 1 };
36: PetscViewer view;
37: PetscReal norm;
38: PetscBool verbose,FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;
40: SlepcInitialize(&argc,&argv,(char*)0,help);
41: PetscOptionsGetInt(NULL,"-n",&n,NULL);
42: PetscOptionsGetInt(NULL,"-kx",&kx,NULL);
43: PetscOptionsGetInt(NULL,"-lx",&lx,NULL);
44: PetscOptionsGetInt(NULL,"-ky",&ky,NULL);
45: PetscOptionsGetInt(NULL,"-ly",&ly,NULL);
46: PetscOptionsHasName(NULL,"-verbose",&verbose);
47: PetscPrintf(PETSC_COMM_WORLD,"Test BV projection (n=%D).\n",n);
48: PetscPrintf(PETSC_COMM_WORLD,"X has %D active columns (%D leading columns).\n",kx,lx);
49: PetscPrintf(PETSC_COMM_WORLD,"Y has %D active columns (%D leading columns).\n",ky,ly);
51: /* Set up viewer */
52: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view);
53: if (verbose) {
54: PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB);
55: }
57: /* Create non-symmetric matrix G (Toeplitz) */
58: MatCreate(PETSC_COMM_WORLD,&G);
59: MatSetSizes(G,PETSC_DECIDE,PETSC_DECIDE,n,n);
60: MatSetFromOptions(G);
61: MatSetUp(G);
62: PetscObjectSetName((PetscObject)G,"G");
64: MatGetOwnershipRange(G,&Istart,&Iend);
65: for (i=Istart;i<Iend;i++) {
66: col[0]=i-1; col[1]=i; col[2]=i+1; col[3]=i+2; col[4]=i+3;
67: if (i==0) {
68: MatSetValues(G,1,&i,PetscMin(4,n-i),col+1,value+1,INSERT_VALUES);
69: } else {
70: MatSetValues(G,1,&i,PetscMin(5,n-i+1),col,value,INSERT_VALUES);
71: }
72: }
73: MatAssemblyBegin(G,MAT_FINAL_ASSEMBLY);
74: MatAssemblyEnd(G,MAT_FINAL_ASSEMBLY);
75: if (verbose) {
76: MatView(G,view);
77: }
79: /* Create symmetric matrix B (1-D Laplacian) */
80: MatCreate(PETSC_COMM_WORLD,&B);
81: MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n);
82: MatSetFromOptions(B);
83: MatSetUp(B);
84: PetscObjectSetName((PetscObject)B,"B");
86: MatGetOwnershipRange(B,&Istart,&Iend);
87: if (Istart==0) FirstBlock=PETSC_TRUE;
88: if (Iend==n) LastBlock=PETSC_TRUE;
89: value[0]=-1.0; value[1]=2.0; value[2]=-1.0;
90: for (i=(FirstBlock? Istart+1: Istart); i<(LastBlock? Iend-1: Iend); i++) {
91: col[0]=i-1; col[1]=i; col[2]=i+1;
92: MatSetValues(B,1,&i,3,col,value,INSERT_VALUES);
93: }
94: if (LastBlock) {
95: i=n-1; col[0]=n-2; col[1]=n-1;
96: MatSetValues(B,1,&i,2,col,value,INSERT_VALUES);
97: }
98: if (FirstBlock) {
99: i=0; col[0]=0; col[1]=1; value[0]=2.0; value[1]=-1.0;
100: MatSetValues(B,1,&i,2,col,value,INSERT_VALUES);
101: }
102: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
103: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
104: MatGetVecs(B,&t,NULL);
105: if (verbose) {
106: MatView(B,view);
107: }
109: /* Create BV object X */
110: BVCreate(PETSC_COMM_WORLD,&X);
111: PetscObjectSetName((PetscObject)X,"X");
112: BVSetSizesFromVec(X,t,kx+2); /* two extra columns to test active columns */
113: BVSetFromOptions(X);
115: /* Fill X entries */
116: for (j=0;j<kx+2;j++) {
117: BVGetColumn(X,j,&v);
118: VecZeroEntries(v);
119: for (i=0;i<4;i++) {
120: if (i+j<n) {
121: VecSetValue(v,i+j,(PetscScalar)(3*i+j-2),INSERT_VALUES);
122: }
123: }
124: VecAssemblyBegin(v);
125: VecAssemblyEnd(v);
126: BVRestoreColumn(X,j,&v);
127: }
128: if (verbose) {
129: BVView(X,view);
130: }
132: /* Duplicate BV object and store Z=G*X */
133: BVDuplicate(X,&Z);
134: PetscObjectSetName((PetscObject)Z,"Z");
135: BVSetActiveColumns(X,0,kx);
136: BVMatMult(X,G,Z);
137: BVSetActiveColumns(X,lx,kx);
138: BVSetActiveColumns(Z,lx,kx);
140: /* Create BV object Y */
141: BVCreate(PETSC_COMM_WORLD,&Y);
142: PetscObjectSetName((PetscObject)Y,"Y");
143: BVSetSizesFromVec(Y,t,ky+1);
144: BVSetFromOptions(Y);
145: BVSetActiveColumns(Y,ly,ky);
147: /* Fill Y entries */
148: for (j=0;j<ky+1;j++) {
149: BVGetColumn(Y,j,&v);
150: VecSet(v,(PetscScalar)(j+1)/4.0);
151: BVRestoreColumn(Y,j,&v);
152: }
153: if (verbose) {
154: BVView(Y,view);
155: }
157: /* Test BVMatProject for non-symmetric matrix G */
158: MatCreateSeqDense(PETSC_COMM_SELF,ky,kx,NULL,&H0);
159: PetscObjectSetName((PetscObject)H0,"H0");
160: BVMatProject(X,G,Y,H0);
161: if (verbose) {
162: MatView(H0,view);
163: }
165: /* Test BVMatProject with previously stored G*X */
166: MatCreateSeqDense(PETSC_COMM_SELF,ky,kx,NULL,&H1);
167: PetscObjectSetName((PetscObject)H1,"H1");
168: BVMatProject(Z,NULL,Y,H1);
169: if (verbose) {
170: MatView(H1,view);
171: }
173: /* Check that H0 and H1 are equal */
174: MatAXPY(H0,-1.0,H1,SAME_NONZERO_PATTERN);
175: MatNorm(H0,NORM_1,&norm);
176: if (norm<10*PETSC_MACHINE_EPSILON) {
177: PetscPrintf(PETSC_COMM_WORLD,"||H0-H1|| < 10*eps\n");
178: } else {
179: PetscPrintf(PETSC_COMM_WORLD,"||H0-H1||=%g\n",(double)norm);
180: }
181: MatDestroy(&H0);
182: MatDestroy(&H1);
184: /* Test BVMatProject for symmetric matrix B with orthogonal projection */
185: MatCreateSeqDense(PETSC_COMM_SELF,kx,kx,NULL,&H0);
186: PetscObjectSetName((PetscObject)H0,"H0");
187: BVMatProject(X,B,X,H0);
188: if (verbose) {
189: MatView(H0,view);
190: }
192: /* Repeat previous test with symmetry flag set */
193: MatSetOption(B,MAT_HERMITIAN,PETSC_TRUE);
194: MatCreateSeqDense(PETSC_COMM_SELF,kx,kx,NULL,&H1);
195: PetscObjectSetName((PetscObject)H1,"H1");
196: BVMatProject(X,B,X,H1);
197: if (verbose) {
198: MatView(H1,view);
199: }
201: /* Check that H0 and H1 are equal */
202: MatAXPY(H0,-1.0,H1,SAME_NONZERO_PATTERN);
203: MatNorm(H0,NORM_1,&norm);
204: if (norm<10*PETSC_MACHINE_EPSILON) {
205: PetscPrintf(PETSC_COMM_WORLD,"||H0-H1|| < 10*eps\n");
206: } else {
207: PetscPrintf(PETSC_COMM_WORLD,"||H0-H1||=%g\n",(double)norm);
208: }
209: MatDestroy(&H0);
210: MatDestroy(&H1);
212: BVDestroy(&X);
213: BVDestroy(&Y);
214: BVDestroy(&Z);
215: MatDestroy(&B);
216: MatDestroy(&G);
217: VecDestroy(&t);
218: SlepcFinalize();
219: return 0;
220: }