Actual source code: test5.c

slepc-3.18.2 2023-01-26
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;

 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 with indefinite inner product (n=%" PetscInt_FMT ", k=%" PetscInt_FMT ").\n",n,k);

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

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

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

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

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

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

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

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

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

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

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

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

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

137: /*TEST

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

157: TEST*/