Actual source code: test12.c

slepc-3.11.2 2019-07-30
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2019, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */

 11: static char help[] = "Test block orthogonalization on a rank-deficient BV.\n\n";

 13: #include <slepcbv.h>

 15: int main(int argc,char **argv)
 16: {
 18:   BV             X,Z;
 19:   Mat            M,R;
 20:   Vec            v,w,t;
 21:   PetscInt       i,j,n=20,k=8;
 22:   PetscViewer    view;
 23:   PetscBool      verbose;
 24:   PetscReal      norm;
 25:   PetscScalar    alpha;

 27:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 28:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 29:   PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL);
 30:   PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
 31:   PetscPrintf(PETSC_COMM_WORLD,"Test BV block orthogonalization (length %D, k=%D).\n",n,k);
 32:   if (k<6) SETERRQ(PETSC_COMM_WORLD,1,"k must be at least 6");

 34:   /* Create template vector */
 35:   VecCreate(PETSC_COMM_WORLD,&t);
 36:   VecSetSizes(t,PETSC_DECIDE,n);
 37:   VecSetFromOptions(t);

 39:   /* Create BV object X */
 40:   BVCreate(PETSC_COMM_WORLD,&X);
 41:   PetscObjectSetName((PetscObject)X,"X");
 42:   BVSetSizesFromVec(X,t,k);
 43:   BVSetFromOptions(X);

 45:   /* Set up viewer */
 46:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view);
 47:   if (verbose) {
 48:     PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB);
 49:   }

 51:   /* Fill X entries (first half) */
 52:   for (j=0;j<k/2;j++) {
 53:     BVGetColumn(X,j,&v);
 54:     VecSet(v,0.0);
 55:     for (i=0;i<=n/2;i++) {
 56:       if (i+j<n) {
 57:         alpha = (3.0*i+j-2)/(2*(i+j+1));
 58:         VecSetValue(v,i+j,alpha,INSERT_VALUES);
 59:       }
 60:     }
 61:     VecAssemblyBegin(v);
 62:     VecAssemblyEnd(v);
 63:     BVRestoreColumn(X,j,&v);
 64:   }

 66:   /* make middle column linearly dependent wrt columns 0 and 1 */
 67:   BVCopyColumn(X,0,j);
 68:   BVGetColumn(X,j,&v);
 69:   BVGetColumn(X,1,&w);
 70:   VecAXPY(v,0.5,w);
 71:   BVRestoreColumn(X,1,&w);
 72:   BVRestoreColumn(X,j,&v);
 73:   j++;

 75:   /* Fill X entries (second half) */
 76:   for (;j<k-1;j++) {
 77:     BVGetColumn(X,j,&v);
 78:     VecSet(v,0.0);
 79:     for (i=0;i<=n/2;i++) {
 80:       if (i+j<n) {
 81:         alpha = (3.0*i+j-2)/(2*(i+j+1));
 82:         VecSetValue(v,i+j,alpha,INSERT_VALUES);
 83:       }
 84:     }
 85:     VecAssemblyBegin(v);
 86:     VecAssemblyEnd(v);
 87:     BVRestoreColumn(X,j,&v);
 88:   }

 90:   /* make middle column linearly dependent wrt columns 1 and k/2+1 */
 91:   BVCopyColumn(X,1,j);
 92:   BVGetColumn(X,j,&v);
 93:   BVGetColumn(X,k/2+1,&w);
 94:   VecAXPY(v,-1.2,w);
 95:   BVRestoreColumn(X,k/2+1,&w);
 96:   BVRestoreColumn(X,j,&v);

 98:   if (verbose) {
 99:     BVView(X,view);
100:   }

102:   /* Create a copy on Z */
103:   BVDuplicate(X,&Z);
104:   PetscObjectSetName((PetscObject)Z,"Z");
105:   BVCopy(X,Z);

107:   /* Test BVOrthogonalize */
108:   MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&R);
109:   PetscObjectSetName((PetscObject)R,"R");
110:   BVOrthogonalize(X,R);
111:   if (verbose) {
112:     BVView(X,view);
113:     MatView(R,view);
114:   }

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

128:   /* Check residual */
129:   BVMult(Z,-1.0,1.0,X,R);
130:   BVNorm(Z,NORM_FROBENIUS,&norm);
131:   if (norm<100*PETSC_MACHINE_EPSILON) {
132:     PetscPrintf(PETSC_COMM_WORLD,"Residual ||X-QR|| < 100*eps\n");
133:   } else {
134:     PetscPrintf(PETSC_COMM_WORLD,"Residual ||X-QR||: %g\n",(double)norm);
135:   }

137:   MatDestroy(&R);
138:   MatDestroy(&M);
139:   BVDestroy(&X);
140:   BVDestroy(&Z);
141:   VecDestroy(&t);
142:   SlepcFinalize();
143:   return ierr;
144: }

146: /*TEST

148:    test:
149:       suffix: 1
150:       nsize: 1
151:       args: -bv_orthog_block gs -bv_type {{vecs contiguous svec mat}shared output}
152:       output_file: output/test12_1.out

154: TEST*/