Actual source code: test2.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 SVD with different builds with a matrix loaded from a file"
 12:   " (matrices available in PETSc's distribution).\n\n";

 14: #include <slepcsvd.h>

 16: int main(int argc,char **argv)
 17: {
 18:   Mat            A;               /* operator matrix */
 19:   SVD            svd;             /* singular value problem solver context */
 20:   char           filename[PETSC_MAX_PATH_LEN];
 21:   const char     *prefix,*scalar,*ints,*floats;
 22:   PetscReal      tol=1000*PETSC_MACHINE_EPSILON;
 23:   PetscViewer    viewer;

 26:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;

 28:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 29:         Load the matrix for which the SVD must be computed
 30:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 31: #if defined(PETSC_USE_COMPLEX)
 32:   prefix = "nh";
 33:   scalar = "complex";
 34: #else
 35:   prefix = "ns";
 36:   scalar = "real";
 37: #endif
 38: #if defined(PETSC_USE_64BIT_INDICES)
 39:   ints   = "int64";
 40: #else
 41:   ints   = "int32";
 42: #endif
 43: #if defined(PETSC_USE_REAL_DOUBLE)
 44:   floats = "float64";
 45: #elif defined(PETSC_USE_REAL_SINGLE)
 46:   floats = "float32";
 47: #endif

 49:   PetscSNPrintf(filename,PETSC_MAX_PATH_LEN,"%s/share/petsc/datafiles/matrices/%s-%s-%s-%s",PETSC_DIR,prefix,scalar,ints,floats);
 50:   PetscPrintf(PETSC_COMM_WORLD,"\nReading matrix from binary file...\n\n");
 51:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);
 52:   MatCreate(PETSC_COMM_WORLD,&A);
 53:   MatSetFromOptions(A);
 54:   MatLoad(A,viewer);
 55:   PetscViewerDestroy(&viewer);

 57:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 58:                      Create the SVD solver
 59:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 60:   SVDCreate(PETSC_COMM_WORLD,&svd);
 61:   SVDSetOperator(svd,A);
 62:   SVDSetTolerances(svd,tol,PETSC_DEFAULT);
 63:   SVDSetFromOptions(svd);

 65:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 66:                 Compute the singular triplets and display solution
 67:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 68:   SVDSolve(svd);
 69:   SVDErrorView(svd,SVD_ERROR_RELATIVE,NULL);
 70:   SVDDestroy(&svd);
 71:   MatDestroy(&A);
 72:   SlepcFinalize();
 73:   return ierr;
 74: }

 76: /*TEST

 78:    build:
 79:       requires: !__float128

 81:    testset:
 82:       args: -svd_nsv 7
 83:       test:
 84:          args: -svd_type {{lanczos trlanczos cross cyclic lapack}}
 85:          requires: double
 86:       test:
 87:          suffix: 1_primme
 88:          args: -svd_type primme
 89:          requires: primme
 90:          output_file: output/test2_1.out

 92: TEST*/