Actual source code: test5.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 view and monitor functionality.\n\n";

 13: #include <slepcsvd.h>

 15: int main(int argc,char **argv)
 16: {
 17:   Mat            A;
 18:   SVD            svd;
 19:   PetscInt       n=6,Istart,Iend,i;

 22:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 23:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 24:   PetscPrintf(PETSC_COMM_WORLD,"\nSVD of diagonal matrix, n=%D\n\n",n);

 26:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 27:         Generate the matrix
 28:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 29:   MatCreate(PETSC_COMM_WORLD,&A);
 30:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
 31:   MatSetFromOptions(A);
 32:   MatSetUp(A);
 33:   MatGetOwnershipRange(A,&Istart,&Iend);
 34:   for (i=Istart;i<Iend;i++) {
 35:     MatSetValue(A,i,i,i+1,INSERT_VALUES);
 36:   }
 37:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 38:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 40:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 41:                      Create the SVD solver
 42:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 43:   SVDCreate(PETSC_COMM_WORLD,&svd);
 44:   PetscObjectSetName((PetscObject)svd,"svd");
 45:   SVDSetOperator(svd,A);
 46:   SVDSetFromOptions(svd);

 48:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 49:                 Compute the singular triplets and display solution
 50:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 51:   SVDSolve(svd);
 52:   SVDErrorView(svd,SVD_ERROR_RELATIVE,NULL);
 53:   SVDDestroy(&svd);
 54:   MatDestroy(&A);
 55:   SlepcFinalize();
 56:   return ierr;
 57: }

 59: /*TEST

 61:    test:
 62:       suffix: 1
 63:       args: -svd_error_relative ::ascii_info_detail -svd_view_values -svd_monitor_conv -svd_error_absolute ::ascii_matlab -svd_monitor_all -svd_converged_reason -svd_view
 64:       filter: grep -v "tolerance" | grep -v "problem type" | sed -e "s/1.999999/2.000000/" | sed -e "s/[0-9]\.[0-9]*e-\([0-9]*\)/removed/g"

 66: TEST*/