Actual source code: test8.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 orthogonalization with selected columns.\n\n";

 24: #include <slepcbv.h>

 28: int main(int argc,char **argv)
 29: {
 31:   BV             X;
 32:   Vec            v,t,z;
 33:   PetscInt       i,j,n=20,k=8;
 34:   PetscViewer    view;
 35:   PetscBool      verbose,*which;
 36:   PetscReal      norm;
 37:   PetscScalar    alpha,*pz;

 39:   SlepcInitialize(&argc,&argv,(char*)0,help);
 40:   PetscOptionsGetInt(NULL,"-n",&n,NULL);
 41:   PetscOptionsGetInt(NULL,"-k",&k,NULL);
 42:   PetscOptionsHasName(NULL,"-verbose",&verbose);
 43:   PetscPrintf(PETSC_COMM_WORLD,"Test BV orthogonalization with selected columns of length %D.\n",n);

 45:   /* Create template vector */
 46:   VecCreate(PETSC_COMM_WORLD,&t);
 47:   VecSetSizes(t,PETSC_DECIDE,n);
 48:   VecSetFromOptions(t);

 50:   /* Create BV object X */
 51:   BVCreate(PETSC_COMM_WORLD,&X);
 52:   PetscObjectSetName((PetscObject)X,"X");
 53:   BVSetSizesFromVec(X,t,k);
 54:   BVSetOrthogonalization(X,BV_ORTHOG_MGS,BV_ORTHOG_REFINE_IFNEEDED,PETSC_DEFAULT,BV_ORTHOG_BLOCK_GS);
 55:   BVSetFromOptions(X);

 57:   /* Set up viewer */
 58:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view);
 59:   if (verbose) {
 60:     PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB);
 61:   }

 63:   /* Fill X entries */
 64:   for (j=0;j<k;j++) {
 65:     BVGetColumn(X,j,&v);
 66:     VecSet(v,0.0);
 67:     for (i=0;i<=n/2;i++) {
 68:       if (i+j<n) {
 69:         alpha = (3.0*i+j-2)/(2*(i+j+1));
 70:         VecSetValue(v,i+j,alpha,INSERT_VALUES);
 71:       }
 72:     }
 73:     VecAssemblyBegin(v);
 74:     VecAssemblyEnd(v);
 75:     BVRestoreColumn(X,j,&v);
 76:   }
 77:   if (verbose) {
 78:     BVView(X,view);
 79:   }

 81:   /* Orthonormalize first k-1 columns */
 82:   for (j=0;j<k-1;j++) {
 83:     BVOrthogonalizeColumn(X,j,NULL,&norm,NULL);
 84:     alpha = 1.0/norm;
 85:     BVScaleColumn(X,j,alpha);
 86:   }
 87:   if (verbose) {
 88:     BVView(X,view);
 89:   }

 91:   /* Select odd columns and orthogonalize last column against those only */
 92:   PetscMalloc1(k,&which);
 93:   for (i=0;i<k;i++) which[i] = (i%2)? PETSC_TRUE: PETSC_FALSE;
 94:   BVOrthogonalizeSomeColumn(X,k-1,which,NULL,NULL,NULL);
 95:   PetscFree(which);
 96:   if (verbose) {
 97:     BVView(X,view);
 98:   }

100:   /* Check orthogonality */
101:   PetscPrintf(PETSC_COMM_WORLD,"Orthogonalization coefficients:\n");
102:   VecCreateSeq(PETSC_COMM_SELF,k-1,&z);
103:   PetscObjectSetName((PetscObject)z,"z");
104:   VecGetArray(z,&pz);
105:   BVDotColumn(X,k-1,pz);
106:   for (i=0;i<k-1;i++) {
107:     if (PetscAbsScalar(pz[i])<PETSC_MACHINE_EPSILON) pz[i]=0.0;
108:   }
109:   VecRestoreArray(z,&pz);
110:   VecView(z,view);
111:   VecDestroy(&z);

113:   BVDestroy(&X);
114:   VecDestroy(&t);
115:   SlepcFinalize();
116:   return 0;
117: }