Actual source code: test12.c

slepc-3.14.2 2021-02-01
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2020, 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[] = "SVD problem with user-defined stopping test.\n\n"
 12:   "The command line options are:\n"
 13:   "  -m <m>, where <m> = matrix rows.\n"
 14:   "  -n <n>, where <n> = matrix columns (defaults to m+2).\n\n";

 16: #include <slepcsvd.h>
 17: #include <petsctime.h>

 19: /*
 20:    This example computes the singular values of a rectangular bidiagonal matrix

 22:               |  1  2                     |
 23:               |     1  2                  |
 24:               |        1  2               |
 25:           A = |          .  .             |
 26:               |             .  .          |
 27:               |                1  2       |
 28:               |                   1  2    |
 29:  */

 31: /*
 32:     Function for user-defined stopping test.

 34:     Checks that the computing time has not exceeded the deadline.
 35: */
 36: PetscErrorCode MyStoppingTest(SVD svd,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,SVDConvergedReason *reason,void *ctx)
 37: {
 39:   PetscLogDouble now,deadline = *(PetscLogDouble*)ctx;

 42:   /* check if usual termination conditions are met */
 43:   SVDStoppingBasic(svd,its,max_it,nconv,nev,reason,NULL);
 44:   if (*reason==SVD_CONVERGED_ITERATING) {
 45:     /* check if deadline has expired */
 46:     PetscTime(&now);
 47:     if (now>deadline) *reason = SVD_CONVERGED_USER;
 48:   }
 49:   return(0);
 50: }

 52: int main(int argc,char **argv)
 53: {
 54:   Mat                A;
 55:   SVD                svd;
 56:   SVDConvergedReason reason;
 57:   PetscInt           m=20,n,Istart,Iend,i,col[2],nconv;
 58:   PetscReal          seconds=2.5;     /* maximum time allowed for computation */
 59:   PetscLogDouble     deadline;        /* time to abort computation */
 60:   PetscScalar        value[] = { 1, 2 };
 61:   PetscBool          terse,flg;
 62:   PetscViewer        viewer;
 63:   PetscErrorCode     ierr;

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

 67:   PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
 68:   PetscOptionsGetInt(NULL,NULL,"-n",&n,&flg);
 69:   if (!flg) n=m+2;
 70:   PetscPrintf(PETSC_COMM_WORLD,"\nRectangular bidiagonal matrix, m=%D n=%D\n\n",m,n);
 71:   PetscOptionsGetReal(NULL,NULL,"-seconds",&seconds,NULL);
 72:   deadline = seconds;
 73:   PetscTimeAdd(&deadline);

 75:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 76:         Generate the matrix
 77:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 79:   MatCreate(PETSC_COMM_WORLD,&A);
 80:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,n);
 81:   MatSetFromOptions(A);
 82:   MatSetUp(A);
 83:   MatGetOwnershipRange(A,&Istart,&Iend);
 84:   for (i=Istart;i<Iend;i++) {
 85:     col[0]=i; col[1]=i+1;
 86:     if (i<n-1) {
 87:       MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
 88:     } else if (i==n-1) {
 89:       MatSetValue(A,i,col[0],value[0],INSERT_VALUES);
 90:     }
 91:   }
 92:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 93:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 95:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 96:         Compute singular values
 97:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 99:   SVDCreate(PETSC_COMM_WORLD,&svd);
100:   SVDSetOperator(svd,A);
101:   SVDSetWhichSingularTriplets(svd,SVD_SMALLEST);
102:   SVDSetTolerances(svd,PETSC_DEFAULT,1000);
103:   SVDSetType(svd,SVDTRLANCZOS);
104:   SVDSetStoppingTestFunction(svd,MyStoppingTest,&deadline,NULL);
105:   SVDSetFromOptions(svd);

107:   /* call the solver */
108:   SVDSolve(svd);

110:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
111:                     Display solution and clean up
112:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

114:   /* show detailed info unless -terse option is given by user */
115:   PetscOptionsHasName(NULL,NULL,"-terse",&terse);
116:   if (terse) {
117:     SVDErrorView(svd,SVD_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
118:   } else {
119:     PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
120:     PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
121:     SVDGetConvergedReason(svd,&reason);
122:     if (reason!=SVD_CONVERGED_USER) {
123:       SVDConvergedReasonView(svd,viewer);
124:       SVDErrorView(svd,SVD_ERROR_RELATIVE,viewer);
125:     } else {
126:       SVDGetConverged(svd,&nconv);
127:       PetscViewerASCIIPrintf(viewer,"SVD solve finished with %D converged eigenpairs; reason=%s\n",nconv,SVDConvergedReasons[reason]);
128:     }
129:     PetscViewerPopFormat(viewer);
130:   }
131:   SVDDestroy(&svd);
132:   MatDestroy(&A);
133:   SlepcFinalize();
134:   return ierr;
135: }

137: /*TEST

139:    test:
140:       suffix: 1
141:       args: -m 750 -seconds 0.1 -svd_max_it 10000

143: TEST*/