Actual source code: test3.c

slepc-3.6.3 2016-03-29
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2015, 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;
 35:   PetscScalar    alpha;
 36:   PetscReal      nrm;
 37:   PetscViewer    view;
 38:   PetscBool      verbose;

 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:   for (i=Istart;i<Iend;i++) {
 55:     if (i>0) { MatSetValue(B,i,i-1,-1.0,INSERT_VALUES); }
 56:     if (i<n-1) { MatSetValue(B,i,i+1,-1.0,INSERT_VALUES); }
 57:     MatSetValue(B,i,i,2.0,INSERT_VALUES);
 58:   }
 59:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 60:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
 61:   MatCreateVecs(B,&t,NULL);

 63:   /* Create BV object X */
 64:   BVCreate(PETSC_COMM_WORLD,&X);
 65:   PetscObjectSetName((PetscObject)X,"X");
 66:   BVSetSizesFromVec(X,t,k);
 67:   BVSetFromOptions(X);
 68:   BVSetMatrix(X,B,PETSC_FALSE);

 70:   /* Set up viewer */
 71:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view);
 72:   if (verbose) {
 73:     PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB);
 74:   }

 76:   /* Fill X entries */
 77:   for (j=0;j<k;j++) {
 78:     BVGetColumn(X,j,&v);
 79:     VecSet(v,0.0);
 80:     for (i=0;i<4;i++) {
 81:       if (i+j<n) {
 82:         VecSetValue(v,i+j,(PetscScalar)(3*i+j-2),INSERT_VALUES);
 83:       }
 84:     }
 85:     VecAssemblyBegin(v);
 86:     VecAssemblyEnd(v);
 87:     BVRestoreColumn(X,j,&v);
 88:   }
 89:   if (verbose) {
 90:     MatView(B,view);
 91:     BVView(X,view);
 92:   }

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

 98:   /* Test BVOrthogonalizeColumn */
 99:   for (j=0;j<k;j++) {
100:     BVOrthogonalizeColumn(X,j,NULL,&nrm,NULL);
101:     alpha = 1.0/nrm;
102:     BVScaleColumn(X,j,alpha);
103:   }
104:   if (verbose) {
105:     BVView(X,view);
106:   }

108:   /* Check orthogonality */
109:   MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&M);
110:   BVDot(X,X,M);
111:   MatShift(M,-1.0);
112:   MatNorm(M,NORM_1,&nrm);
113:   if (nrm<100*PETSC_MACHINE_EPSILON) {
114:     PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality < 100*eps\n");
115:   } else {
116:     PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g\n",(double)nrm);
117:   }

119:   BVDestroy(&X);
120:   MatDestroy(&M);
121:   MatDestroy(&B);
122:   VecDestroy(&t);
123:   SlepcFinalize();
124:   return 0;
125: }