Actual source code: test17.c

slepc-3.17.2 2022-08-09
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 bi-orthogonalization functions.\n\n";

 13: #include <slepcbv.h>

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

 26:   SlepcInitialize(&argc,&argv,(char*)0,help);
 27:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 28:   PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL);
 29:   PetscOptionsGetReal(NULL,NULL,"-condn",&condn,NULL);
 31:   PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
 32:   PetscPrintf(PETSC_COMM_WORLD,"Test BV bi-orthogonalization with %" PetscInt_FMT " columns of length %" PetscInt_FMT ".\n",k,n);
 33:   if (condn>1.0) PetscPrintf(PETSC_COMM_WORLD," - Using random BVs with condition number = %g\n",(double)condn);

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

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

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

 50:   /* Fill X entries */
 51:   if (condn==1.0) {
 52:     for (j=0;j<k;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: #if defined(PETSC_USE_COMPLEX)
 59:           alpha += (1.2*i+j-2)/(0.1*(i+j+1))*PETSC_i;
 60: #endif
 61:           VecSetValue(v,i+j,alpha,INSERT_VALUES);
 62:         }
 63:       }
 64:       VecAssemblyBegin(v);
 65:       VecAssemblyEnd(v);
 66:       BVRestoreColumn(X,j,&v);
 67:     }
 68:   } else BVSetRandomCond(X,condn);
 69:   if (verbose) BVView(X,view);

 71:   /* Create Y and fill its entries */
 72:   BVDuplicate(X,&Y);
 73:   PetscObjectSetName((PetscObject)Y,"Y");
 74:   if (condn==1.0) {
 75:     for (j=0;j<k;j++) {
 76:       BVGetColumn(Y,j,&v);
 77:       VecSet(v,0.0);
 78:       for (i=PetscMin(n,2+(2*j)%6);i<PetscMin(n,6+(3*j)%9);i++) {
 79:         if (i%5 && i!=j) {
 80:           alpha = (1.5*i+j)/(2.2*(i-j));
 81:           VecSetValue(v,i+j,alpha,INSERT_VALUES);
 82:         }
 83:       }
 84:       VecAssemblyBegin(v);
 85:       VecAssemblyEnd(v);
 86:       BVRestoreColumn(Y,j,&v);
 87:     }
 88:   } else BVSetRandomCond(Y,condn);
 89:   if (verbose) BVView(Y,view);

 91:   /* Test BVBiorthonormalizeColumn */
 92:   for (j=0;j<k;j++) BVBiorthonormalizeColumn(X,Y,j,NULL);
 93:   if (verbose) {
 94:     BVView(X,view);
 95:     BVView(Y,view);
 96:   }

 98:   /* Check orthogonality */
 99:   MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&M);
100:   PetscObjectSetName((PetscObject)M,"M");
101:   BVDot(X,Y,M);
102:   if (verbose) MatView(M,view);
103:   MatShift(M,-1.0);
104:   MatNorm(M,NORM_1,&norm);
105:   if (norm<200*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"Level of bi-orthogonality < 200*eps\n");
106:   else PetscPrintf(PETSC_COMM_WORLD,"Level of bi-orthogonality: %g\n",(double)norm);

108:   MatDestroy(&M);
109:   BVDestroy(&X);
110:   BVDestroy(&Y);
111:   VecDestroy(&t);
112:   SlepcFinalize();
113:   return 0;
114: }

116: /*TEST

118:    testset:
119:       output_file: output/test17_1.out
120:       test:
121:          suffix: 1
122:          args: -bv_type {{vecs contiguous svec mat}} -bv_orthog_type cgs
123:          requires: double
124:       test:
125:          suffix: 1_cuda
126:          args: -bv_type svec -vec_type cuda -bv_orthog_type cgs
127:          requires: cuda
128:       test:
129:          suffix: 2
130:          args: -bv_type {{vecs contiguous svec mat}} -bv_orthog_type mgs
131:       test:
132:          suffix: 2_cuda
133:          args: -bv_type svec -vec_type cuda -bv_orthog_type mgs
134:          requires: cuda

136: TEST*/