Actual source code: test15.c

slepc-3.16.3 2022-04-11
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, 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[] = "Tests user interface for TRLANCZOS with GSVD.\n\n"
 12:   "The command line options are:\n"
 13:   "  -m <m>, where <m> = number of rows of A.\n"
 14:   "  -n <n>, where <n> = number of columns of A.\n"
 15:   "  -p <p>, where <p> = number of rows of B.\n\n";

 17: #include <slepcsvd.h>

 19: int main(int argc,char **argv)
 20: {
 21:   Mat                 A,B;
 22:   SVD                 svd;
 23:   KSP                 ksp;
 24:   PC                  pc;
 25:   PetscInt            m=15,n=20,p=21,i,j,d,Istart,Iend;
 26:   PetscReal           keep;
 27:   PetscBool           flg,lock;
 28:   SVDTRLanczosGBidiag bidiag;
 29:   PetscErrorCode      ierr;

 31:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 32:   PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
 33:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 34:   PetscOptionsGetInt(NULL,NULL,"-p",&p,NULL);
 35:   PetscPrintf(PETSC_COMM_WORLD,"\nGeneralized singular value decomposition, (%D+%D)x%D\n\n",m,p,n);

 37:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 38:                      Generate the matrices
 39:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 41:   MatCreate(PETSC_COMM_WORLD,&A);
 42:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,n);
 43:   MatSetFromOptions(A);
 44:   MatSetUp(A);

 46:   MatGetOwnershipRange(A,&Istart,&Iend);
 47:   for (i=Istart;i<Iend;i++) {
 48:     if (i>0 && i-1<n) { MatSetValue(A,i,i-1,-1.0,INSERT_VALUES); }
 49:     if (i+1<n) { MatSetValue(A,i,i+1,-1.0,INSERT_VALUES); }
 50:     if (i<n) { MatSetValue(A,i,i,2.0,INSERT_VALUES); }
 51:     if (i>n) { MatSetValue(A,i,n-1,1.0,INSERT_VALUES); }
 52:   }
 53:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 54:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 56:   MatCreate(PETSC_COMM_WORLD,&B);
 57:   MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,p,n);
 58:   MatSetFromOptions(B);
 59:   MatSetUp(B);

 61:   MatGetOwnershipRange(B,&Istart,&Iend);
 62:   d = PetscMax(0,n-p);
 63:   for (i=Istart;i<Iend;i++) {
 64:     for (j=0;j<=PetscMin(i,n-1);j++) {
 65:       MatSetValue(B,i,j+d,1.0,INSERT_VALUES);
 66:     }
 67:   }
 68:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 69:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);

 71:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 72:          Create the singular value solver, set options and solve
 73:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 75:   SVDCreate(PETSC_COMM_WORLD,&svd);
 76:   SVDSetOperators(svd,A,B);
 77:   SVDSetDimensions(svd,4,PETSC_DEFAULT,PETSC_DEFAULT);
 78:   SVDSetConvergenceTest(svd,SVD_CONV_NORM);

 80:   SVDSetType(svd,SVDTRLANCZOS);
 81:   SVDTRLanczosSetGBidiag(svd,SVD_TRLANCZOS_GBIDIAG_UPPER);

 83:   /* create a standalone KSP with appropriate settings */
 84:   KSPCreate(PETSC_COMM_WORLD,&ksp);
 85:   KSPSetType(ksp,KSPLSQR);
 86:   KSPGetPC(ksp,&pc);
 87:   PCSetType(pc,PCNONE);
 88:   KSPSetTolerances(ksp,PETSC_SMALL,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
 89:   KSPSetFromOptions(ksp);
 90:   SVDTRLanczosSetKSP(svd,ksp);
 91:   SVDTRLanczosSetRestart(svd,0.4);
 92:   SVDTRLanczosSetLocking(svd,PETSC_TRUE);

 94:   SVDSetFromOptions(svd);

 96:   PetscObjectTypeCompare((PetscObject)svd,SVDTRLANCZOS,&flg);
 97:   if (flg) {
 98:     SVDTRLanczosGetGBidiag(svd,&bidiag);
 99:     PetscPrintf(PETSC_COMM_WORLD,"TRLANCZOS: using %s bidiagonalization\n",SVDTRLanczosGBidiags[bidiag]);
100:     SVDTRLanczosGetRestart(svd,&keep);
101:     SVDTRLanczosGetLocking(svd,&lock);
102:     PetscPrintf(PETSC_COMM_WORLD,"TRLANCZOS: restarting parameter %.2f %s\n",(double)keep,lock?"(locking)":"");
103:   }

105:   SVDSolve(svd);
106:   SVDErrorView(svd,SVD_ERROR_RELATIVE,NULL);

108:   /* Free work space */
109:   SVDDestroy(&svd);
110:   KSPDestroy(&ksp);
111:   MatDestroy(&A);
112:   MatDestroy(&B);
113:   SlepcFinalize();
114:   return ierr;
115: }

117: /*TEST

119:    test:
120:       suffix: 1
121:       requires: !single

123:    test:
124:       suffix: 2
125:       args: -m 6 -n 12 -p 12 -svd_trlanczos_gbidiag {{single upper lower}}
126:       filter: grep -v "TRLANCZOS: using"
127:       requires: !single

129: TEST*/