Actual source code: test12.c

slepc-3.18.3 2023-03-24
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, 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: {
 17:   BV             X,Z;
 18:   Mat            M,R;
 19:   Vec            v,w,t;
 20:   PetscInt       i,j,n=20,k=8;
 21:   PetscViewer    view;
 22:   PetscBool      verbose;
 23:   PetscReal      norm;
 24:   PetscScalar    alpha;

 27:   SlepcInitialize(&argc,&argv,(char*)0,help);
 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 %" PetscInt_FMT ", k=%" PetscInt_FMT ").\n",n,k);

 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) PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB);

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

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

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

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

 96:   if (verbose) BVView(X,view);

 98:   /* Create a copy on Z */
 99:   BVDuplicate(X,&Z);
100:   PetscObjectSetName((PetscObject)Z,"Z");
101:   BVCopy(X,Z);

103:   /* Test BVOrthogonalize */
104:   MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&R);
105:   PetscObjectSetName((PetscObject)R,"R");
106:   BVOrthogonalize(X,R);
107:   if (verbose) {
108:     BVView(X,view);
109:     MatView(R,view);
110:   }

112:   /* Check orthogonality */
113:   MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&M);
114:   MatShift(M,1.0);   /* set leading part to identity */
115:   BVDot(X,X,M);
116:   MatShift(M,-1.0);
117:   MatNorm(M,NORM_1,&norm);
118:   if (norm<100*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality < 100*eps\n");
119:   else PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g\n",(double)norm);

121:   /* Check residual */
122:   BVMult(Z,-1.0,1.0,X,R);
123:   BVNorm(Z,NORM_FROBENIUS,&norm);
124:   if (norm<100*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"Residual ||X-QR|| < 100*eps\n");
125:   else PetscPrintf(PETSC_COMM_WORLD,"Residual ||X-QR||: %g\n",(double)norm);

127:   MatDestroy(&R);
128:   MatDestroy(&M);
129:   BVDestroy(&X);
130:   BVDestroy(&Z);
131:   VecDestroy(&t);
132:   SlepcFinalize();
133:   return 0;
134: }

136: /*TEST

138:    test:
139:       suffix: 1
140:       nsize: 1
141:       args: -bv_orthog_block gs -bv_type {{vecs contiguous svec mat}shared output}
142:       output_file: output/test12_1.out

144: TEST*/