Actual source code: test2.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 BV orthogonalization functions.\n\n";

 13: #include <slepcbv.h>

 15: int main(int argc,char **argv)
 16: {
 18:   BV             X,Y,Z;
 19:   Mat            M,R;
 20:   Vec            v,t,e;
 21:   PetscInt       i,j,n=20,k=8;
 22:   PetscViewer    view;
 23:   PetscBool      verbose;
 24:   PetscReal      norm,condn=1.0;
 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:   PetscOptionsGetReal(NULL,NULL,"-condn",&condn,NULL);
 31:   if (condn<1.0) SETERRQ(PETSC_COMM_WORLD,1,"The condition number must be > 1");
 32:   PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
 33:   PetscPrintf(PETSC_COMM_WORLD,"Test BV orthogonalization with %D columns of length %D.\n",k,n);
 34:   if (condn>1.0) {
 35:     PetscPrintf(PETSC_COMM_WORLD," - Using a random BV with condition number = %g\n",(double)condn);
 36:   }

 38:   /* Create template vector */
 39:   VecCreate(PETSC_COMM_WORLD,&t);
 40:   VecSetSizes(t,PETSC_DECIDE,n);
 41:   VecSetFromOptions(t);

 43:   /* Create BV object X */
 44:   BVCreate(PETSC_COMM_WORLD,&X);
 45:   PetscObjectSetName((PetscObject)X,"X");
 46:   BVSetSizesFromVec(X,t,k);
 47:   BVSetFromOptions(X);

 49:   /* Set up viewer */
 50:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view);
 51:   if (verbose) {
 52:     PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB);
 53:   }

 55:   /* Fill X entries */
 56:   if (condn==1.0) {
 57:     for (j=0;j<k;j++) {
 58:       BVGetColumn(X,j,&v);
 59:       VecSet(v,0.0);
 60:       for (i=0;i<=n/2;i++) {
 61:         if (i+j<n) {
 62:           alpha = (3.0*i+j-2)/(2*(i+j+1));
 63:           VecSetValue(v,i+j,alpha,INSERT_VALUES);
 64:         }
 65:       }
 66:       VecAssemblyBegin(v);
 67:       VecAssemblyEnd(v);
 68:       BVRestoreColumn(X,j,&v);
 69:     }
 70:   } else {
 71:     BVSetRandomCond(X,condn);
 72:   }
 73:   if (verbose) {
 74:     BVView(X,view);
 75:   }

 77:   /* Create copies on Y and Z */
 78:   BVDuplicate(X,&Y);
 79:   PetscObjectSetName((PetscObject)Y,"Y");
 80:   BVCopy(X,Y);
 81:   BVDuplicate(X,&Z);
 82:   PetscObjectSetName((PetscObject)Z,"Z");
 83:   BVCopy(X,Z);

 85:   /* Test BVOrthogonalizeColumn */
 86:   for (j=0;j<k;j++) {
 87:     BVOrthogonalizeColumn(X,j,NULL,&norm,NULL);
 88:     alpha = 1.0/norm;
 89:     BVScaleColumn(X,j,alpha);
 90:   }
 91:   if (verbose) {
 92:     BVView(X,view);
 93:   }

 95:   /* Check orthogonality */
 96:   MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&M);
 97:   BVDot(X,X,M);
 98:   MatShift(M,-1.0);
 99:   MatNorm(M,NORM_1,&norm);
100:   if (norm<100*PETSC_MACHINE_EPSILON) {
101:     PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality < 100*eps\n");
102:   } else {
103:     PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g\n",(double)norm);
104:   }

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

115:   /* Check orthogonality */
116:   BVDot(Y,Y,M);
117:   MatShift(M,-1.0);
118:   MatNorm(M,NORM_1,&norm);
119:   if (norm<100*PETSC_MACHINE_EPSILON) {
120:     PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality < 100*eps\n");
121:   } else {
122:     PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g\n",(double)norm);
123:   }

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

134:   /* Test BVOrthogonalizeVec */
135:   VecDuplicate(t,&e);
136:   VecSet(e,1.0);
137:   BVOrthogonalizeVec(X,e,NULL,&norm,NULL);
138:   PetscPrintf(PETSC_COMM_WORLD,"Norm of ones(n,1) after orthogonalizing against X: %g\n",(double)norm);

140:   MatDestroy(&M);
141:   MatDestroy(&R);
142:   BVDestroy(&X);
143:   BVDestroy(&Y);
144:   BVDestroy(&Z);
145:   VecDestroy(&e);
146:   VecDestroy(&t);
147:   SlepcFinalize();
148:   return ierr;
149: }

151: /*TEST

153:    test:
154:       suffix: 1
155:       nsize: 1
156:       args: -bv_type {{vecs contiguous svec mat}shared output} -bv_orthog_type cgs
157:       output_file: output/test2_1.out

159:    test:
160:       suffix: 1_cuda
161:       nsize: 1
162:       args: -bv_type svec -vec_type cuda -bv_orthog_type cgs
163:       requires: cuda
164:       output_file: output/test2_1.out

166:    test:
167:       suffix: 2
168:       nsize: 1
169:       args: -bv_type {{vecs contiguous svec mat}shared output} -bv_orthog_type mgs
170:       output_file: output/test2_1.out

172:    test:
173:       suffix: 2_cuda
174:       nsize: 1
175:       args: -bv_type svec -vec_type cuda -bv_orthog_type mgs
176:       requires: cuda
177:       output_file: output/test2_1.out

179:    test:
180:       suffix: 3
181:       nsize: 1
182:       args: -bv_type {{vecs contiguous svec mat}shared output} -condn 1e8
183:       requires: !single
184:       filter: grep -v "against"
185:       output_file: output/test2_3.out

187: TEST*/