Actual source code: test14.c

slepc-3.11.2 2019-07-30
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2019, 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: {
 18:   BV             X;
 19:   Mat            A,B,M;
 20:   PetscInt       i,j,n=20,k=8,Istart,Iend;
 21:   PetscViewer    view;
 22:   PetscBool      verbose;
 23:   PetscReal      norm;
 24:   PetscScalar    alpha;

 26:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 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 created from a dense Mat (length %D, k=%D).\n",n,k);

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

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

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

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

 71:   /* Test BVOrthogonalize */
 72:   BVOrthogonalize(X,NULL);
 73:   if (verbose) {
 74:     BVView(X,view);
 75:   }

 77:   /* Check orthogonality */
 78:   MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&M);
 79:   MatShift(M,1.0);   /* set leading part to identity */
 80:   BVDot(X,X,M);
 81:   MatShift(M,-1.0);
 82:   MatNorm(M,NORM_1,&norm);
 83:   if (norm<100*PETSC_MACHINE_EPSILON) {
 84:     PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality < 100*eps\n");
 85:   } else {
 86:     PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g\n",(double)norm);
 87:   }

 89:   MatDestroy(&M);
 90:   MatDestroy(&A);
 91:   MatDestroy(&B);
 92:   BVDestroy(&X);
 93:   SlepcFinalize();
 94:   return ierr;
 95: }

 97: /*TEST

 99:    test:
100:       suffix: 1
101:       nsize: 2
102:       args: -bv_type {{vecs contiguous svec mat}shared output}
103:       output_file: output/test14_1.out

105: TEST*/