Actual source code: test2.c

slepc-3.22.1 2024-10-28
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 BV orthogonalization functions.\n\n";

 13: #include <slepcbv.h>

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

 26:   PetscFunctionBeginUser;
 27:   PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
 28:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 29:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL));
 30:   PetscCall(PetscOptionsGetReal(NULL,NULL,"-condn",&condn,NULL));
 31:   PetscCheck(condn>=1.0,PETSC_COMM_WORLD,PETSC_ERR_USER_INPUT,"The condition number must be > 1");
 32:   PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose));
 33:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test BV orthogonalization with %" PetscInt_FMT " columns of length %" PetscInt_FMT ".\n",k,n));
 34:   if (condn>1.0) PetscCall(PetscPrintf(PETSC_COMM_WORLD," - Using a random BV with condition number = %g\n",(double)condn));

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

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

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

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

 69:   /* Create copies on Y and Z */
 70:   PetscCall(BVDuplicate(X,&Y));
 71:   PetscCall(PetscObjectSetName((PetscObject)Y,"Y"));
 72:   PetscCall(BVCopy(X,Y));
 73:   PetscCall(BVDuplicate(X,&Z));
 74:   PetscCall(PetscObjectSetName((PetscObject)Z,"Z"));
 75:   PetscCall(BVCopy(X,Z));

 77:   /* Test BVOrthogonalizeColumn */
 78:   for (j=0;j<k;j++) {
 79:     PetscCall(BVOrthogonalizeColumn(X,j,NULL,&norm,NULL));
 80:     alpha = 1.0/norm;
 81:     PetscCall(BVScaleColumn(X,j,alpha));
 82:   }
 83:   if (verbose) PetscCall(BVView(X,view));

 85:   /* Check orthogonality */
 86:   PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&M));
 87:   PetscCall(BVDot(X,X,M));
 88:   PetscCall(MatShift(M,-1.0));
 89:   PetscCall(MatNorm(M,NORM_1,&norm));
 90:   if (norm<100*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality < 100*eps\n"));
 91:   else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g\n",(double)norm));

 93:   /* Test BVOrthogonalize */
 94:   PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&R));
 95:   PetscCall(PetscObjectSetName((PetscObject)R,"R"));
 96:   PetscCall(BVOrthogonalize(Y,R));
 97:   if (verbose) {
 98:     PetscCall(BVView(Y,view));
 99:     PetscCall(MatView(R,view));
100:   }

102:   /* Check orthogonality */
103:   PetscCall(BVDot(Y,Y,M));
104:   PetscCall(MatShift(M,-1.0));
105:   PetscCall(MatNorm(M,NORM_1,&norm));
106:   if (norm<100*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality < 100*eps\n"));
107:   else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g\n",(double)norm));

109:   /* Check residual */
110:   PetscCall(BVMult(Z,-1.0,1.0,Y,R));
111:   PetscCall(BVNorm(Z,NORM_FROBENIUS,&norm));
112:   if (norm<100*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Residual ||X-QR|| < 100*eps\n"));
113:   else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Residual ||X-QR||: %g\n",(double)norm));

115:   /* Test BVOrthogonalizeVec */
116:   PetscCall(VecDuplicate(t,&e));
117:   PetscCall(VecSet(e,1.0));
118:   PetscCall(BVOrthogonalizeVec(X,e,NULL,&norm,NULL));
119:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Norm of ones(n,1) after orthogonalizing against X: %g\n",(double)norm));

121:   PetscCall(MatDestroy(&M));
122:   PetscCall(MatDestroy(&R));
123:   PetscCall(BVDestroy(&X));
124:   PetscCall(BVDestroy(&Y));
125:   PetscCall(BVDestroy(&Z));
126:   PetscCall(VecDestroy(&e));
127:   PetscCall(VecDestroy(&t));
128:   PetscCall(SlepcFinalize());
129:   return 0;
130: }

132: /*TEST

134:    testset:
135:       output_file: output/test2_1.out
136:       test:
137:          suffix: 1
138:          args: -bv_type {{vecs contiguous svec mat}shared output} -bv_orthog_type cgs
139:       test:
140:          suffix: 1_cuda
141:          args: -bv_type {{svec mat}} -vec_type cuda -bv_orthog_type cgs
142:          requires: cuda
143:       test:
144:          suffix: 1_hip
145:          args: -bv_type {{svec mat}} -vec_type hip -bv_orthog_type cgs
146:          requires: hip
147:       test:
148:          suffix: 2
149:          args: -bv_type {{vecs contiguous svec mat}shared output} -bv_orthog_type mgs
150:       test:
151:          suffix: 2_cuda
152:          args: -bv_type {{svec mat}} -vec_type cuda -bv_orthog_type mgs
153:          requires: cuda
154:       test:
155:          suffix: 2_hip
156:          args: -bv_type {{svec mat}} -vec_type hip -bv_orthog_type mgs
157:          requires: hip

159:    test:
160:       suffix: 3
161:       nsize: 1
162:       args: -bv_type {{vecs contiguous svec mat}shared output} -condn 1e8
163:       requires: !single
164:       filter: grep -v "against"
165:       output_file: output/test2_3.out

167: TEST*/