Actual source code: test4.c

slepc-3.7.4 2017-05-17
Report Typos and Errors
  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, changing the number of active columns.\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,kx=6,lx=3,ky=5,ly=2;
 35:   PetscScalar    *q,*z;
 36:   PetscReal      nrm;
 37:   PetscViewer    view;
 38:   PetscBool      verbose,trans;

 40:   SlepcInitialize(&argc,&argv,(char*)0,help);
 41:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 42:   PetscOptionsGetInt(NULL,NULL,"-kx",&kx,NULL);
 43:   PetscOptionsGetInt(NULL,NULL,"-lx",&lx,NULL);
 44:   PetscOptionsGetInt(NULL,NULL,"-ky",&ky,NULL);
 45:   PetscOptionsGetInt(NULL,NULL,"-ly",&ly,NULL);
 46:   PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
 47:   PetscPrintf(PETSC_COMM_WORLD,"First BV with %D active columns (%D leading columns) of dimension %D.\n",kx,lx,n);
 48:   PetscPrintf(PETSC_COMM_WORLD,"Second BV with %D active columns (%D leading columns) of dimension %D.\n",ky,ly,n);

 50:   /* Create template vector */
 51:   VecCreate(PETSC_COMM_WORLD,&t);
 52:   VecSetSizes(t,PETSC_DECIDE,n);
 53:   VecSetFromOptions(t);

 55:   /* Create BV object X */
 56:   BVCreate(PETSC_COMM_WORLD,&X);
 57:   PetscObjectSetName((PetscObject)X,"X");
 58:   BVSetSizesFromVec(X,t,kx+2);  /* two extra columns to test active columns */
 59:   BVSetFromOptions(X);
 60:   BVSetActiveColumns(X,lx,kx);

 62:   /* Set up viewer */
 63:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view);
 64:   if (verbose) {
 65:     PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB);
 66:   }

 68:   /* Fill X entries */
 69:   for (j=0;j<kx+2;j++) {
 70:     BVGetColumn(X,j,&v);
 71:     VecSet(v,0.0);
 72:     for (i=0;i<4;i++) {
 73:       if (i+j<n) {
 74:         VecSetValue(v,i+j,(PetscScalar)(3*i+j-2),INSERT_VALUES);
 75:       }
 76:     }
 77:     VecAssemblyBegin(v);
 78:     VecAssemblyEnd(v);
 79:     BVRestoreColumn(X,j,&v);
 80:   }
 81:   if (verbose) {
 82:     BVView(X,view);
 83:   }

 85:   /* Create BV object Y */
 86:   BVCreate(PETSC_COMM_WORLD,&Y);
 87:   PetscObjectSetName((PetscObject)Y,"Y");
 88:   BVSetSizesFromVec(Y,t,ky+1);
 89:   BVSetFromOptions(Y);
 90:   BVSetActiveColumns(Y,ly,ky);

 92:   /* Fill Y entries */
 93:   for (j=0;j<ky+1;j++) {
 94:     BVGetColumn(Y,j,&v);
 95:     VecSet(v,(PetscScalar)(j+1)/4.0);
 96:     BVRestoreColumn(Y,j,&v);
 97:   }
 98:   if (verbose) {
 99:     BVView(Y,view);
100:   }

102:   /* Create Mat */
103:   MatCreateSeqDense(PETSC_COMM_SELF,kx,ky,NULL,&Q);
104:   PetscObjectSetName((PetscObject)Q,"Q");
105:   MatDenseGetArray(Q,&q);
106:   for (i=0;i<kx;i++)
107:     for (j=0;j<ky;j++)
108:       q[i+j*kx] = (i<j)? 2.0: -0.5;
109:   MatDenseRestoreArray(Q,&q);
110:   if (verbose) {
111:     MatView(Q,NULL);
112:   }

114:   /* Test BVMult */
115:   BVMult(Y,2.0,1.0,X,Q);
116:   if (verbose) {
117:     PetscPrintf(PETSC_COMM_WORLD,"After BVMult - - - - - - - - -\n");
118:     BVView(Y,view);
119:   }

121:   /* Test BVMultVec */
122:   BVGetColumn(Y,0,&v);
123:   PetscMalloc1(kx-lx,&z);
124:   z[0] = 2.0;
125:   for (i=1;i<kx-lx;i++) z[i] = -0.5*z[i-1];
126:   BVMultVec(X,-1.0,1.0,v,z);
127:   PetscFree(z);
128:   BVRestoreColumn(Y,0,&v);
129:   if (verbose) {
130:     PetscPrintf(PETSC_COMM_WORLD,"After BVMultVec - - - - - - -\n");
131:     BVView(Y,view);
132:   }

134:   /* Test BVDot */
135:   MatCreateSeqDense(PETSC_COMM_SELF,ky,kx,NULL,&M);
136:   PetscObjectSetName((PetscObject)M,"M");
137:   BVDot(X,Y,M);
138:   if (verbose) {
139:     PetscPrintf(PETSC_COMM_WORLD,"After BVDot - - - - - - - - -\n");
140:     MatView(M,NULL);
141:   }

143:   /* Test BVDotVec */
144:   BVGetColumn(Y,0,&v);
145:   PetscMalloc1(kx-lx,&z);
146:   BVDotVec(X,v,z);
147:   BVRestoreColumn(Y,0,&v);
148:   if (verbose) {
149:     PetscPrintf(PETSC_COMM_WORLD,"After BVDotVec - - - - - - -\n");
150:     VecCreateSeqWithArray(PETSC_COMM_SELF,1,kx-lx,z,&v);
151:     PetscObjectSetName((PetscObject)v,"z");
152:     VecView(v,view);
153:     VecDestroy(&v);
154:   }
155:   PetscFree(z);

157:   /* Test BVMultInPlace and BVScale */
158:   PetscOptionsHasName(NULL,NULL,"-trans",&trans);
159:   if (trans) {
160:     Mat Qt;
161:     MatTranspose(Q,MAT_INITIAL_MATRIX,&Qt);
162:     BVMultInPlaceTranspose(X,Qt,lx+1,ky);
163:     MatDestroy(&Qt);
164:   } else {
165:     BVMultInPlace(X,Q,lx+1,ky);
166:   }
167:   BVScale(X,2.0);
168:   if (verbose) {
169:     PetscPrintf(PETSC_COMM_WORLD,"After BVMultInPlace - - - - -\n");
170:     BVView(X,view);
171:   }

173:   /* Test BVNorm */
174:   BVNormColumn(X,lx,NORM_2,&nrm);
175:   PetscPrintf(PETSC_COMM_WORLD,"2-Norm or X[%D] = %g\n",lx,(double)nrm);
176:   BVNorm(X,NORM_FROBENIUS,&nrm);
177:   PetscPrintf(PETSC_COMM_WORLD,"Frobenius Norm or X = %g\n",(double)nrm);

179:   BVDestroy(&X);
180:   BVDestroy(&Y);
181:   MatDestroy(&Q);
182:   MatDestroy(&M);
183:   VecDestroy(&t);
184:   SlepcFinalize();
185:   return ierr;
186: }