Actual source code: test2.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 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=PETSC_SMALL;
 23:   PetscViewer    viewer;

 25:   SlepcInitialize(&argc,&argv,(char*)0,help);

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

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

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

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

 75: /*TEST

 77:    build:
 78:       requires: !__float128

 80:    test:
 81:       args: -svd_nsv 7 -svd_type {{lanczos trlanczos cross cyclic lapack randomized}}
 82:       requires: !single

 84:    testset:
 85:       args: -svd_nsv 7 -svd_mpd 11 -svd_type primme
 86:       requires: primme !single
 87:       output_file: output/test2_1.out
 88:       test:
 89:          suffix: 1_primme
 90:       test:
 91:          suffix: 1_primme_args
 92:          args: -svd_primme_blocksize 2 -svd_primme_method hybrid

 94: TEST*/