Actual source code: test5.c

slepc-3.17.0 2022-03-31
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 operations with indefinite inner product.\n\n";

 13: #include <slepcbv.h>

 15: int main(int argc,char **argv)
 16: {
 17:   Vec            t,v,w,omega;
 18:   Mat            B,M;
 19:   BV             X,Y;
 20:   PetscInt       i,j,n=10,k=5,l,Istart,Iend;
 21:   PetscScalar    alpha;
 22:   PetscReal      nrm;
 23:   PetscViewer    view;
 24:   PetscBool      verbose;

 26:   SlepcInitialize(&argc,&argv,(char*)0,help);
 27:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 28:   PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL);
 29:   PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
 30:   PetscPrintf(PETSC_COMM_WORLD,"Test BV with indefinite inner product (n=%" PetscInt_FMT ", k=%" PetscInt_FMT ").\n",n,k);

 32:   /* Create inner product matrix (standard involutionary permutation) */
 33:   MatCreate(PETSC_COMM_WORLD,&B);
 34:   MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n);
 35:   MatSetFromOptions(B);
 36:   MatSetUp(B);
 37:   PetscObjectSetName((PetscObject)B,"B");

 39:   MatGetOwnershipRange(B,&Istart,&Iend);
 40:   for (i=Istart;i<Iend;i++) MatSetValue(B,i,n-i-1,1.0,INSERT_VALUES);
 41:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 42:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
 43:   MatCreateVecs(B,&t,NULL);

 45:   /* Create BV object X */
 46:   BVCreate(PETSC_COMM_WORLD,&X);
 47:   PetscObjectSetName((PetscObject)X,"X");
 48:   BVSetSizesFromVec(X,t,k);
 49:   BVSetFromOptions(X);
 50:   BVSetMatrix(X,B,PETSC_TRUE);

 52:   /* Set up viewer */
 53:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view);
 54:   if (verbose) PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB);

 56:   /* Fill X entries */
 57:   l = -3;
 58:   for (j=0;j<k;j++) {
 59:     BVGetColumn(X,j,&v);
 60:     VecSet(v,-1.0);
 61:     for (i=0;i<n/2;i++) {
 62:       if (i+j<n) {
 63:         l = (l + 3*i+j-2) % n;
 64:         VecSetValue(v,i+j,(PetscScalar)l,INSERT_VALUES);
 65:       }
 66:     }
 67:     VecAssemblyBegin(v);
 68:     VecAssemblyEnd(v);
 69:     BVRestoreColumn(X,j,&v);
 70:   }
 71:   if (verbose) {
 72:     MatView(B,view);
 73:     BVView(X,view);
 74:   }

 76:   /* Test BVNormColumn */
 77:   BVNormColumn(X,0,NORM_2,&nrm);
 78:   PetscPrintf(PETSC_COMM_WORLD,"B-Norm of X[0] = %g\n",(double)nrm);

 80:   /* Test BVOrthogonalizeColumn */
 81:   for (j=0;j<k;j++) {
 82:     BVOrthogonalizeColumn(X,j,NULL,&nrm,NULL);
 83:     alpha = 1.0/nrm;
 84:     BVScaleColumn(X,j,alpha);
 85:   }
 86:   if (verbose) BVView(X,view);

 88:   /* Create a copy on Y */
 89:   BVDuplicate(X,&Y);
 90:   PetscObjectSetName((PetscObject)Y,"Y");
 91:   BVCopy(X,Y);

 93:   /* Check orthogonality */
 94:   MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&M);
 95:   BVDot(Y,Y,M);
 96:   VecCreateSeq(PETSC_COMM_SELF,k,&omega);
 97:   BVGetSignature(Y,omega);
 98:   VecScale(omega,-1.0);
 99:   MatDiagonalSet(M,omega,ADD_VALUES);
100:   MatNorm(M,NORM_1,&nrm);
101:   if (nrm<100*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality < 100*eps\n");
102:   else PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g\n",(double)nrm);

104:   /* Test BVSetSignature */
105:   VecScale(omega,-1.0);
106:   BVSetSignature(Y,omega);
107:   VecDestroy(&omega);

109:   /* Test BVApplyMatrix */
110:   VecDuplicate(t,&w);
111:   BVGetColumn(X,0,&v);
112:   BVApplyMatrix(X,v,w);
113:   BVApplyMatrix(X,w,t);
114:   VecAXPY(t,-1.0,v);
115:   BVRestoreColumn(X,0,&v);
116:   VecNorm(t,NORM_2,&nrm);

119:   BVApplyMatrixBV(X,Y);
120:   BVGetColumn(Y,0,&v);
121:   VecAXPY(w,-1.0,v);
122:   BVRestoreColumn(Y,0,&v);
123:   VecNorm(w,NORM_2,&nrm);

126:   BVDestroy(&X);
127:   BVDestroy(&Y);
128:   MatDestroy(&M);
129:   MatDestroy(&B);
130:   VecDestroy(&w);
131:   VecDestroy(&t);
132:   SlepcFinalize();
133:   return 0;
134: }

136: /*TEST

138:    testset:
139:       output_file: output/test5_1.out
140:       args: -bv_orthog_refine always
141:       test:
142:          suffix: 1
143:          args: -bv_type {{vecs contiguous svec mat}shared output}
144:       test:
145:          suffix: 1_cuda
146:          args: -bv_type svec -mat_type aijcusparse
147:          requires: cuda
148:       test:
149:          suffix: 2
150:          args: -bv_type {{vecs contiguous svec mat}shared output} -bv_orthog_type mgs
151:       test:
152:          suffix: 2_cuda
153:          args: -bv_type svec -mat_type aijcusparse -bv_orthog_type mgs
154:          requires: cuda

156: TEST*/