Actual source code: test14.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 created from a dense Mat.\n\n";

 13: #include <slepcbv.h>

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

 25:   SlepcInitialize(&argc,&argv,(char*)0,help);
 26:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 27:   PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL);
 28:   PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
 29:   PetscPrintf(PETSC_COMM_WORLD,"Test BV created from a dense Mat (length %" PetscInt_FMT ", k=%" PetscInt_FMT ").\n",n,k);

 31:   /* Create dense matrix */
 32:   MatCreate(PETSC_COMM_WORLD,&A);
 33:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,k);
 34:   MatSetType(A,MATDENSE);
 35:   MatSetUp(A);
 36:   MatGetOwnershipRange(A,&Istart,&Iend);
 37:   for (j=0;j<k;j++) {
 38:     for (i=0;i<=n/2;i++) {
 39:       if (i+j<n && i>=Istart && i<Iend) {
 40:         alpha = (3.0*i+j-2)/(2*(i+j+1));
 41:         MatSetValue(A,i+j,j,alpha,INSERT_VALUES);
 42:       }
 43:     }
 44:   }
 45:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 46:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 48:   /* Create BV object X */
 49:   BVCreateFromMat(A,&X);
 50:   BVSetFromOptions(X);
 51:   PetscObjectSetName((PetscObject)X,"X");

 53:   /* Set up viewer */
 54:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&view);
 55:   if (verbose) {
 56:     PetscViewerPushFormat(view,PETSC_VIEWER_ASCII_MATLAB);
 57:     BVView(X,view);
 58:   }

 60:   /* Test BVCreateMat */
 61:   BVCreateMat(X,&B);
 62:   MatAXPY(B,-1.0,A,SAME_NONZERO_PATTERN);
 63:   MatNorm(B,NORM_1,&norm);
 64:   if (norm<100*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"Norm of difference < 100*eps\n");
 65:   else PetscPrintf(PETSC_COMM_WORLD,"Norm of difference: %g\n",(double)norm);

 67:   /* Test BVOrthogonalize */
 68:   BVOrthogonalize(X,NULL);
 69:   if (verbose) BVView(X,view);

 71:   /* Check orthogonality */
 72:   MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&M);
 73:   MatShift(M,1.0);   /* set leading part to identity */
 74:   BVDot(X,X,M);
 75:   MatShift(M,-1.0);
 76:   MatNorm(M,NORM_1,&norm);
 77:   if (norm<100*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality < 100*eps\n");
 78:   else PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g\n",(double)norm);

 80:   MatDestroy(&M);
 81:   MatDestroy(&A);
 82:   MatDestroy(&B);
 83:   BVDestroy(&X);
 84:   SlepcFinalize();
 85:   return 0;
 86: }

 88: /*TEST

 90:    test:
 91:       suffix: 1
 92:       nsize: 2
 93:       args: -bv_type {{vecs contiguous svec mat}shared output}
 94:       output_file: output/test14_1.out

 96: TEST*/