Actual source code: test10.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 split reductions in BV.\n\n";
24: #include <slepcbv.h>
28: int main(int argc,char **argv)
29: {
31: Vec t,v,w,y,z,zsplit;
32: BV X;
33: PetscInt i,j,n=10,k=5;
34: PetscScalar *zarray;
35: PetscReal nrm;
37: SlepcInitialize(&argc,&argv,(char*)0,help);
38: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
39: PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL);
40: if (k<3) SETERRQ(PETSC_COMM_SELF,1,"Should specify at least k=3 columns");
41: PetscPrintf(PETSC_COMM_WORLD,"BV split ops (%D columns of dimension %D).\n",k,n);
43: /* Create template vector */
44: VecCreate(PETSC_COMM_WORLD,&t);
45: VecSetSizes(t,PETSC_DECIDE,n);
46: VecSetFromOptions(t);
47: VecDuplicate(t,&v);
48: VecSet(v,1.0);
50: /* Create BV object X */
51: BVCreate(PETSC_COMM_WORLD,&X);
52: PetscObjectSetName((PetscObject)X,"X");
53: BVSetSizesFromVec(X,t,k);
54: BVSetFromOptions(X);
56: /* Fill X entries */
57: for (j=0;j<k;j++) {
58: BVGetColumn(X,j,&w);
59: VecSet(w,0.0);
60: for (i=0;i<4;i++) {
61: if (i+j<n) {
62: VecSetValue(w,i+j,(PetscScalar)(3*i+j-2),INSERT_VALUES);
63: }
64: }
65: VecAssemblyBegin(w);
66: VecAssemblyEnd(w);
67: BVRestoreColumn(X,j,&w);
68: }
70: /* Use regular operations */
71: VecCreateSeq(PETSC_COMM_SELF,k+6,&z);
72: VecGetArray(z,&zarray);
73: BVGetColumn(X,0,&w);
74: VecDot(w,v,zarray);
75: BVRestoreColumn(X,0,&w);
76: BVDotVec(X,v,zarray+1);
77: BVDotColumn(X,2,zarray+1+k);
79: BVGetColumn(X,1,&y);
80: VecNorm(y,NORM_2,&nrm);
81: BVRestoreColumn(X,1,&y);
82: zarray[k+3] = nrm;
83: BVNormVec(X,v,NORM_2,&nrm);
84: zarray[k+4] = nrm;
85: BVNormColumn(X,0,NORM_2,&nrm);
86: zarray[k+5] = nrm;
87: VecRestoreArray(z,&zarray);
89: /* Use split operations */
90: VecCreateSeq(PETSC_COMM_SELF,k+6,&zsplit);
91: VecGetArray(zsplit,&zarray);
92: BVGetColumn(X,0,&w);
93: VecDotBegin(w,v,zarray);
94: BVDotVecBegin(X,v,zarray+1);
95: BVDotColumnBegin(X,2,zarray+1+k);
96: VecDotEnd(w,v,zarray);
97: BVRestoreColumn(X,0,&w);
98: BVDotVecEnd(X,v,zarray+1);
99: BVDotColumnEnd(X,2,zarray+1+k);
101: BVGetColumn(X,1,&y);
102: VecNormBegin(y,NORM_2,&nrm);
103: BVNormVecBegin(X,v,NORM_2,&nrm);
104: BVNormColumnBegin(X,0,NORM_2,&nrm);
105: VecNormEnd(y,NORM_2,&nrm);
106: BVRestoreColumn(X,1,&y);
107: zarray[k+3] = nrm;
108: BVNormVecEnd(X,v,NORM_2,&nrm);
109: zarray[k+4] = nrm;
110: BVNormColumnEnd(X,0,NORM_2,&nrm);
111: zarray[k+5] = nrm;
112: VecRestoreArray(zsplit,&zarray);
114: /* Show difference */
115: VecAXPY(z,-1.0,zsplit);
116: VecNorm(z,NORM_1,&nrm);
117: PetscSynchronizedPrintf(PETSC_COMM_WORLD,"%g\n",(double)nrm);
118: PetscSynchronizedFlush(PETSC_COMM_WORLD,NULL);
120: BVDestroy(&X);
121: VecDestroy(&t);
122: VecDestroy(&v);
123: VecDestroy(&z);
124: VecDestroy(&zsplit);
125: SlepcFinalize();
126: return ierr;
127: }