Actual source code: test3.c

slepc-3.5.2 2014-10-10
Report Typos and Errors
  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 operations with non-standard 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:   BV             X;
 34:   PetscInt       i,j,n=10,k=5,Istart,Iend,col[3];
 35:   PetscScalar    value[3],alpha;
 36:   PetscReal      nrm;
 37:   PetscViewer    view;
 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,"-k",&k,NULL);
 43:   PetscOptionsHasName(NULL,"-verbose",&verbose);
 44:   PetscPrintf(PETSC_COMM_WORLD,"Test BV with non-standard inner product (n=%D, k=%D).\n",n,k);

 46:   /* Create inner product matrix */
 47:   MatCreate(PETSC_COMM_WORLD,&B);
 48:   MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n);
 49:   MatSetFromOptions(B);
 50:   MatSetUp(B);
 51:   PetscObjectSetName((PetscObject)B,"B");

 53:   MatGetOwnershipRange(B,&Istart,&Iend);
 54:   if (Istart==0) FirstBlock=PETSC_TRUE;
 55:   if (Iend==n) LastBlock=PETSC_TRUE;
 56:   value[0]=-1.0; value[1]=2.0; value[2]=-1.0;
 57:   for (i=(FirstBlock? Istart+1: Istart); i<(LastBlock? Iend-1: Iend); i++) {
 58:     col[0]=i-1; col[1]=i; col[2]=i+1;
 59:     MatSetValues(B,1,&i,3,col,value,INSERT_VALUES);
 60:   }
 61:   if (LastBlock) {
 62:     i=n-1; col[0]=n-2; col[1]=n-1;
 63:     MatSetValues(B,1,&i,2,col,value,INSERT_VALUES);
 64:   }
 65:   if (FirstBlock) {
 66:     i=0; col[0]=0; col[1]=1; value[0]=2.0; value[1]=-1.0;
 67:     MatSetValues(B,1,&i,2,col,value,INSERT_VALUES);
 68:   }
 69:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 70:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
 71:   MatGetVecs(B,&t,NULL);

 73:   /* Create BV object X */
 74:   BVCreate(PETSC_COMM_WORLD,&X);
 75:   PetscObjectSetName((PetscObject)X,"X");
 76:   BVSetSizesFromVec(X,t,k);
 77:   BVSetFromOptions(X);
 78:   BVSetMatrix(X,B,PETSC_FALSE);

 80:   /* Set up viewer */
 81:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view);
 82:   if (verbose) {
 83:     PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB);
 84:   }

 86:   /* Fill X entries */
 87:   for (j=0;j<k;j++) {
 88:     BVGetColumn(X,j,&v);
 89:     VecZeroEntries(v);
 90:     for (i=0;i<4;i++) {
 91:       if (i+j<n) {
 92:         VecSetValue(v,i+j,(PetscScalar)(3*i+j-2),INSERT_VALUES);
 93:       }
 94:     }
 95:     VecAssemblyBegin(v);
 96:     VecAssemblyEnd(v);
 97:     BVRestoreColumn(X,j,&v);
 98:   }
 99:   if (verbose) {
100:     MatView(B,view);
101:     BVView(X,view);
102:   }

104:   /* Test BVNormColumn */
105:   BVNormColumn(X,0,NORM_2,&nrm);
106:   PetscPrintf(PETSC_COMM_WORLD,"B-Norm or X[0] = %g\n",(double)nrm);

108:   /* Test BVOrthogonalizeColumn */
109:   for (j=0;j<k;j++) {
110:     BVOrthogonalizeColumn(X,j,NULL,&nrm,NULL);
111:     alpha = 1.0/nrm;
112:     BVScaleColumn(X,j,alpha);
113:   }
114:   if (verbose) {
115:     BVView(X,view);
116:   }

118:   /* Check orthogonality */
119:   MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&M);
120:   BVDot(X,X,M);
121:   MatShift(M,-1.0);
122:   MatNorm(M,NORM_1,&nrm);
123:   if (nrm<100*PETSC_MACHINE_EPSILON) {
124:     PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality < 100*eps\n");
125:   } else {
126:     PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g\n",(double)nrm);
127:   }

129:   BVDestroy(&X);
130:   MatDestroy(&M);
131:   MatDestroy(&B);
132:   VecDestroy(&t);
133:   SlepcFinalize();
134:   return 0;
135: }