Actual source code: test1.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 BV operations.\n\n";
24: #include <slepcbv.h>
28: int main(int argc,char **argv)
29: {
31: Vec t,v;
32: Mat Q,M;
33: BV X,Y;
34: PetscInt i,j,n=10,k=5,l=3;
35: PetscScalar *q,*z;
36: PetscReal nrm;
37: PetscViewer view;
38: PetscBool verbose;
40: SlepcInitialize(&argc,&argv,(char*)0,help);
41: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
42: PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL);
43: PetscOptionsGetInt(NULL,NULL,"-l",&l,NULL);
44: PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
45: PetscPrintf(PETSC_COMM_WORLD,"Test BV with %D columns of dimension %D.\n",k,n);
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_INFO_DETAIL);
62: BVView(X,view);
63: PetscViewerPopFormat(view);
64: }
66: /* Fill X entries */
67: for (j=0;j<k;j++) {
68: BVGetColumn(X,j,&v);
69: VecSet(v,0.0);
70: for (i=0;i<4;i++) {
71: if (i+j<n) {
72: VecSetValue(v,i+j,(PetscScalar)(3*i+j-2),INSERT_VALUES);
73: }
74: }
75: VecAssemblyBegin(v);
76: VecAssemblyEnd(v);
77: BVRestoreColumn(X,j,&v);
78: }
79: if (verbose) {
80: BVView(X,view);
81: }
83: /* Create BV object Y */
84: BVCreate(PETSC_COMM_WORLD,&Y);
85: PetscObjectSetName((PetscObject)Y,"Y");
86: BVSetSizesFromVec(Y,t,l);
87: BVSetFromOptions(Y);
89: /* Fill Y entries */
90: for (j=0;j<l;j++) {
91: BVGetColumn(Y,j,&v);
92: VecSet(v,(PetscScalar)(j+1)/4.0);
93: BVRestoreColumn(Y,j,&v);
94: }
95: if (verbose) {
96: BVView(Y,view);
97: }
99: /* Create Mat */
100: MatCreateSeqDense(PETSC_COMM_SELF,k,l,NULL,&Q);
101: PetscObjectSetName((PetscObject)Q,"Q");
102: MatDenseGetArray(Q,&q);
103: for (i=0;i<k;i++)
104: for (j=0;j<l;j++)
105: q[i+j*k] = (i<j)? 2.0: -0.5;
106: MatDenseRestoreArray(Q,&q);
107: if (verbose) {
108: MatView(Q,NULL);
109: }
111: /* Test BVMult */
112: BVMult(Y,2.0,1.0,X,Q);
113: if (verbose) {
114: PetscPrintf(PETSC_COMM_WORLD,"After BVMult - - - - - - - - -\n");
115: BVView(Y,view);
116: }
118: /* Test BVMultVec */
119: BVGetColumn(Y,0,&v);
120: PetscMalloc1(k,&z);
121: z[0] = 2.0;
122: for (i=1;i<k;i++) z[i] = -0.5*z[i-1];
123: BVMultVec(X,-1.0,1.0,v,z);
124: PetscFree(z);
125: BVRestoreColumn(Y,0,&v);
126: if (verbose) {
127: PetscPrintf(PETSC_COMM_WORLD,"After BVMultVec - - - - - - -\n");
128: BVView(Y,view);
129: }
131: /* Test BVDot */
132: MatCreateSeqDense(PETSC_COMM_SELF,l,k,NULL,&M);
133: PetscObjectSetName((PetscObject)M,"M");
134: BVDot(X,Y,M);
135: if (verbose) {
136: PetscPrintf(PETSC_COMM_WORLD,"After BVDot - - - - - - - - -\n");
137: MatView(M,NULL);
138: }
140: /* Test BVDotVec */
141: BVGetColumn(Y,0,&v);
142: PetscMalloc1(k,&z);
143: BVDotVec(X,v,z);
144: BVRestoreColumn(Y,0,&v);
145: if (verbose) {
146: PetscPrintf(PETSC_COMM_WORLD,"After BVDotVec - - - - - - -\n");
147: VecCreateSeqWithArray(PETSC_COMM_SELF,1,k,z,&v);
148: PetscObjectSetName((PetscObject)v,"z");
149: VecView(v,view);
150: VecDestroy(&v);
151: }
152: PetscFree(z);
154: /* Test BVMultInPlace and BVScale */
155: BVMultInPlace(X,Q,1,l);
156: BVScale(X,2.0);
157: if (verbose) {
158: PetscPrintf(PETSC_COMM_WORLD,"After BVMultInPlace - - - - -\n");
159: BVView(X,view);
160: }
162: /* Test BVNorm */
163: BVNormColumn(X,0,NORM_2,&nrm);
164: PetscPrintf(PETSC_COMM_WORLD,"2-Norm or X[0] = %g\n",(double)nrm);
165: BVNorm(X,NORM_FROBENIUS,&nrm);
166: PetscPrintf(PETSC_COMM_WORLD,"Frobenius Norm or X = %g\n",(double)nrm);
168: BVDestroy(&X);
169: BVDestroy(&Y);
170: MatDestroy(&Q);
171: MatDestroy(&M);
172: VecDestroy(&t);
173: SlepcFinalize();
174: return ierr;
175: }