Actual source code: test3.c

slepc-3.19.0 2023-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 non-standard inner product.\n\n";

 13: #include <slepcbv.h>

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

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

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

 40:   PetscCall(MatGetOwnershipRange(B,&Istart,&Iend));
 41:   for (i=Istart;i<Iend;i++) {
 42:     if (i>0) PetscCall(MatSetValue(B,i,i-1,-1.0,INSERT_VALUES));
 43:     if (i<n-1) PetscCall(MatSetValue(B,i,i+1,-1.0,INSERT_VALUES));
 44:     PetscCall(MatSetValue(B,i,i,2.0,INSERT_VALUES));
 45:   }
 46:   PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
 47:   PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
 48:   PetscCall(MatCreateVecs(B,&t,NULL));

 50:   /* Create BV object X */
 51:   PetscCall(BVCreate(PETSC_COMM_WORLD,&X));
 52:   PetscCall(PetscObjectSetName((PetscObject)X,"X"));
 53:   PetscCall(BVSetSizesFromVec(X,t,k));
 54:   PetscCall(BVSetFromOptions(X));
 55:   PetscCall(BVSetMatrix(X,B,PETSC_FALSE));

 57:   /* Set up viewer */
 58:   PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view));
 59:   if (verbose) PetscCall(PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB));

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

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

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

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

 97:   /* Test BVNormVecBegin/End */
 98:   PetscCall(BVGetColumn(X,0,&v));
 99:   PetscCall(BVNormVecBegin(X,v,NORM_1,&nrm));
100:   PetscCall(BVNormVecEnd(X,v,NORM_1,&nrm));
101:   PetscCall(BVRestoreColumn(X,0,&v));
102:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"B-Norm of X[0] = %g\n",(double)nrm));

104:   PetscCall(BVDestroy(&X));
105:   PetscCall(MatDestroy(&M));
106:   PetscCall(MatDestroy(&B));
107:   PetscCall(VecDestroy(&t));
108:   PetscCall(SlepcFinalize());
109:   return 0;
110: }

112: /*TEST

114:    testset:
115:       output_file: output/test3_1.out
116:       test:
117:          suffix: 1
118:          args: -bv_type {{vecs contiguous svec mat}shared output}
119:       test:
120:          suffix: 1_svec_vecs
121:          args: -bv_type svec -bv_matmult vecs
122:       test:
123:          suffix: 1_cuda
124:          args: -bv_type svec -mat_type aijcusparse
125:          requires: cuda
126:       test:
127:          suffix: 2
128:          nsize: 2
129:          args: -bv_type {{vecs contiguous svec mat}shared output}
130:       test:
131:          suffix: 3
132:          nsize: 2
133:          args: -bv_type {{vecs contiguous svec mat}shared output} -bv_orthog_type mgs

135: TEST*/