Actual source code: test6.c

slepc-3.6.1 2015-09-03
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 functions with constraints.\n\n";

 24: #include <slepcbv.h>

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

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

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

 52:   /* Create BV object X */
 53:   BVCreate(PETSC_COMM_WORLD,&X);
 54:   PetscObjectSetName((PetscObject)X,"X");
 55:   BVSetSizesFromVec(X,t,k);
 56:   BVSetFromOptions(X);

 58:   /* Generate constraints and attach them to X */
 59:   if (nc>0) {
 60:     VecDuplicateVecs(t,nc,&C);
 61:     for (j=0;j<nc;j++) {
 62:       for (i=0;i<=j;i++) {
 63:         VecSetValue(C[j],i,1.0,INSERT_VALUES);
 64:       }
 65:       VecAssemblyBegin(C[j]);
 66:       VecAssemblyEnd(C[j]);
 67:     }
 68:     BVInsertConstraints(X,&nc,C);
 69:     VecDestroyVecs(nc,&C);
 70:   }

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

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

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

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

117:   MatDestroy(&M);
118:   BVDestroy(&X);
119:   VecDestroy(&t);
120:   SlepcFinalize();
121:   return 0;
122: }