Actual source code: ex52.c

slepc-3.20.1 2023-11-27
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[] = "Partial hyperbolic singular value decomposition (HSVD) from a file.\n"
 12:   "The command line options are:\n"
 13:   "  -file <filename>, PETSc binary file containing matrix A.\n"
 14:   "  -p <p>, where <p> = number of -1's in signature.\n"
 15:   "  -transpose, to transpose the matrix before doing the computation.\n\n";

 17: #include <slepcsvd.h>

 19: int main(int argc,char **argv)
 20: {
 21:   Mat            A,Omega;         /* operator matrix, signature matrix */
 22:   SVD            svd;             /* singular value problem solver context */
 23:   Mat            At;
 24:   Vec            u,v,vomega,*U,*V;
 25:   MatType        Atype;
 26:   PetscReal      tol,lev1=0.0,lev2=0.0;
 27:   PetscInt       M,N,p=0,i,Istart,Iend,nconv,nsv;
 28:   char           filename[PETSC_MAX_PATH_LEN];
 29:   PetscViewer    viewer;
 30:   PetscBool      flg,terse,skiporth=PETSC_FALSE,transpose=PETSC_FALSE;

 32:   PetscFunctionBeginUser;
 33:   PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));

 35:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 36:         Load matrix that defines the hyperbolic singular value problem
 37:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 39:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nHyperbolic singular value problem stored in file.\n\n"));
 40:   PetscCall(PetscOptionsGetString(NULL,NULL,"-file",filename,sizeof(filename),&flg));
 41:   PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER_INPUT,"Must indicate a file name for matrix A with the -file option");

 43: #if defined(PETSC_USE_COMPLEX)
 44:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Reading COMPLEX matrix from a binary file...\n"));
 45: #else
 46:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Reading REAL matrix from a binary file...\n"));
 47: #endif
 48:   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer));
 49:   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
 50:   PetscCall(MatSetFromOptions(A));
 51:   PetscCall(MatLoad(A,viewer));
 52:   PetscCall(PetscViewerDestroy(&viewer));

 54:   /* transpose the matrix if requested */
 55:   PetscCall(PetscOptionsGetBool(NULL,NULL,"-transpose",&transpose,NULL));
 56:   if (transpose) {
 57:     PetscCall(MatHermitianTranspose(A,MAT_INITIAL_MATRIX,&At));
 58:     PetscCall(MatDestroy(&A));
 59:     A = At;
 60:   }

 62:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 63:                           Create the signature
 64:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 66:   PetscCall(MatGetSize(A,&M,&N));
 67:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-p",&p,&flg));
 68:   PetscCheck(p>=0,PETSC_COMM_WORLD,PETSC_ERR_USER_INPUT,"Parameter p cannot be negative");
 69:   PetscCheck(p<M,PETSC_COMM_WORLD,PETSC_ERR_USER_INPUT,"Parameter p cannot be larger than the number of rows of A");
 70:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n Matrix dimensions: %" PetscInt_FMT "x%" PetscInt_FMT,M,N));
 71:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,", using signature Omega=blkdiag(-I_%" PetscInt_FMT ",I_%" PetscInt_FMT ")\n\n",p,M-p));

 73:   PetscCall(MatCreateVecs(A,NULL,&vomega));
 74:   PetscCall(VecSet(vomega,1.0));
 75:   PetscCall(VecGetOwnershipRange(vomega,&Istart,&Iend));
 76:   for (i=Istart;i<Iend;i++) {
 77:     if (i<p) PetscCall(VecSetValue(vomega,i,-1.0,INSERT_VALUES));
 78:   }
 79:   PetscCall(VecAssemblyBegin(vomega));
 80:   PetscCall(VecAssemblyEnd(vomega));

 82:   PetscCall(MatGetType(A,&Atype));
 83:   PetscCall(MatCreate(PETSC_COMM_WORLD,&Omega));
 84:   PetscCall(MatSetSizes(Omega,PETSC_DECIDE,PETSC_DECIDE,M,M));
 85:   PetscCall(MatSetType(Omega,Atype));
 86:   PetscCall(MatSetUp(Omega));
 87:   PetscCall(MatDiagonalSet(Omega,vomega,INSERT_VALUES));

 89:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 90:           Create the singular value solver and set various options
 91:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 93:   PetscCall(SVDCreate(PETSC_COMM_WORLD,&svd));

 95:   PetscCall(SVDSetOperators(svd,A,NULL));
 96:   PetscCall(SVDSetSignature(svd,vomega));
 97:   PetscCall(SVDSetProblemType(svd,SVD_HYPERBOLIC));

 99:   PetscCall(SVDSetFromOptions(svd));

101:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
102:                  Solve the problem, display solution
103:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

105:   PetscCall(MatCreateVecs(A,&v,&u));
106:   PetscCall(VecSet(u,1.0));
107:   PetscCall(VecSet(v,1.0));
108:   PetscCall(SVDSolve(svd));

110:   /* show detailed info unless -terse option is given by user */
111:   PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
112:   if (terse) PetscCall(SVDErrorView(svd,SVD_ERROR_NORM,NULL));
113:   else {
114:     PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
115:     PetscCall(SVDConvergedReasonView(svd,PETSC_VIEWER_STDOUT_WORLD));
116:     PetscCall(SVDErrorView(svd,SVD_ERROR_NORM,PETSC_VIEWER_STDOUT_WORLD));
117:     PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
118:   }

120:   /* check orthogonality */
121:   PetscCall(PetscOptionsGetBool(NULL,NULL,"-skiporth",&skiporth,NULL));
122:   PetscCall(SVDGetConverged(svd,&nconv));
123:   PetscCall(SVDGetDimensions(svd,&nsv,NULL,NULL));
124:   nconv = PetscMin(nconv,nsv);
125:   if (nconv>0 && !skiporth) {
126:     PetscCall(SVDGetTolerances(svd,&tol,NULL));
127:     PetscCall(VecDuplicateVecs(u,nconv,&U));
128:     PetscCall(VecDuplicateVecs(v,nconv,&V));
129:     for (i=0;i<nconv;i++) PetscCall(SVDGetSingularTriplet(svd,i,NULL,U[i],V[i]));
130:     PetscCall(VecCheckOrthonormality(U,nconv,NULL,nconv,Omega,NULL,&lev1));
131:     PetscCall(VecCheckOrthonormality(V,nconv,NULL,nconv,NULL,NULL,&lev2));
132:     if (lev1+lev2<20*tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality below the tolerance\n"));
133:     else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g (U) %g (V)\n",(double)lev1,(double)lev2));
134:     PetscCall(VecDestroyVecs(nconv,&U));
135:     PetscCall(VecDestroyVecs(nconv,&V));
136:   }
137:   PetscCall(VecDestroy(&u));
138:   PetscCall(VecDestroy(&v));

140:   /* free work space */
141:   PetscCall(SVDDestroy(&svd));
142:   PetscCall(MatDestroy(&A));
143:   PetscCall(MatDestroy(&Omega));
144:   PetscCall(VecDestroy(&vomega));
145:   PetscCall(SlepcFinalize());
146:   return 0;
147: }

149: /*TEST

151:    testset:
152:       args: -file ${DATAFILESPATH}/matrices/real/illc1033.petsc -svd_nsv 62 -p 40 -terse
153:       requires: double !complex datafilespath !defined(PETSC_USE_64BIT_INDICES)
154:       filter: grep -v Reading
155:       output_file: output/ex52_1.out
156:       test:
157:          args: -svd_type cross -svd_cross_explicitmatrix {{0 1}} -svd_implicittranspose {{0 1}}
158:          suffix: 1_cross
159:       test:
160:          args: -svd_type cyclic -svd_cyclic_explicitmatrix {{0 1}} -svd_ncv 300
161:          suffix: 1_cyclic
162:       test:
163:          args: -svd_type trlanczos -svd_trlanczos_explicitmatrix {{0 1}}
164:          suffix: 1_trlanczos
165:       test:
166:          args: -svd_type lapack
167:          suffix: 1_lapack

169:    testset:
170:       args: -file ${DATAFILESPATH}/matrices/real/illc1033.petsc -transpose -svd_nsv 6 -p 130 -terse
171:       requires: double !complex datafilespath !defined(PETSC_USE_64BIT_INDICES)
172:       filter: grep -v Reading
173:       output_file: output/ex52_2.out
174:       test:
175:          args: -svd_type cross -svd_cross_explicitmatrix {{0 1}} -svd_implicittranspose {{0 1}} -svd_ncv 100
176:          suffix: 2_cross
177:       test:
178:          args: -svd_type cyclic -svd_cyclic_explicitmatrix {{0 1}} -svd_ncv 25
179:          suffix: 2_cyclic
180:       test:
181:          args: -svd_type trlanczos -svd_trlanczos_explicitmatrix {{0 1}} -svd_implicittranspose {{0 1}}
182:          suffix: 2_trlanczos

184:    testset:
185:       args: -file ${DATAFILESPATH}/matrices/complex/illc1033.petsc -svd_nsv 62 -p 40 -terse
186:       requires: double complex datafilespath !defined(PETSC_USE_64BIT_INDICES)
187:       filter: grep -v Reading
188:       output_file: output/ex52_1.out
189:       test:
190:          args: -svd_type cross -svd_cross_explicitmatrix {{0 1}}
191:          suffix: 3_cross
192:       test:
193:          args: -svd_type cyclic -svd_cyclic_explicitmatrix {{0 1}} -svd_cyclic_bv_definite_tol 1e-13 -svd_cyclic_st_ksp_type gmres -svd_cyclic_st_pc_type jacobi -svd_ncv 250
194:          suffix: 3_cyclic
195:       test:
196:          args: -svd_type trlanczos -svd_trlanczos_explicitmatrix {{0 1}} -bv_definite_tol 1e-13
197:          suffix: 3_trlanczos
198:       test:
199:          args: -svd_type lapack
200:          suffix: 3_lapack

202:    testset:
203:       args: -file ${DATAFILESPATH}/matrices/complex/illc1033.petsc -transpose -svd_nsv 6 -p 130 -terse
204:       requires: double complex datafilespath !defined(PETSC_USE_64BIT_INDICES)
205:       filter: grep -v Reading
206:       output_file: output/ex52_2.out
207:       test:
208:          args: -svd_type cross -svd_cross_explicitmatrix {{0 1}} -svd_ncv 100 -svd_cross_bv_definite_tol 1e-14
209:          suffix: 4_cross
210:       test:
211:          args: -svd_type cyclic -svd_cyclic_explicitmatrix {{0 1}} -svd_ncv 26
212:          suffix: 4_cyclic
213:       test:
214:          args: -svd_type trlanczos -svd_trlanczos_explicitmatrix {{0 1}}
215:          suffix: 4_trlanczos

217:    testset:
218:       args: -file ${SLEPC_DIR}/share/slepc/datafiles/matrices/rdb200.petsc -svd_smallest -svd_nsv 3 -p 1 -terse
219:       requires: double !complex !defined(PETSC_USE_64BIT_INDICES)
220:       filter: grep -v Reading
221:       output_file: output/ex52_5.out
222:       test:
223:          args: -svd_type cross -svd_max_it 1000 -svd_cross_bv_definite_tol 1e-14
224:          suffix: 5_cross
225:       test:
226:          args: -svd_type cyclic -svd_max_it 4000 -svd_ncv 30 -svd_cyclic_st_ksp_type preonly -svd_cyclic_st_pc_type jacobi -svd_cyclic_bv_definite_tol 1e-14
227:          suffix: 5_cyclic
228:       test:
229:          args: -svd_type trlanczos -svd_max_it 4000 -svd_ncv 28 -bv_definite_tol 1e-14
230:          suffix: 5_trlanczos
231:       test:
232:          args: -svd_type lapack
233:          suffix: 5_lapack

235: TEST*/