Actual source code: test13.c
slepc-3.18.2 2023-01-26
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 interface to external package PRIMME.\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>
18: /*
19: This example computes the singular values of a rectangular bidiagonal matrix
21: | 1 2 |
22: | 1 2 |
23: | 1 2 |
24: A = | . . |
25: | . . |
26: | 1 2 |
27: | 1 2 |
28: */
30: int main(int argc,char **argv)
31: {
32: Mat A;
33: SVD svd;
34: PetscInt m=20,n,Istart,Iend,i,k=6,col[2],bs;
35: PetscScalar value[] = { 1, 2 };
36: PetscBool flg;
37: SVDPRIMMEMethod meth;
40: SlepcInitialize(&argc,&argv,(char*)0,help);
41: PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
42: PetscOptionsGetInt(NULL,NULL,"-n",&n,&flg);
43: if (!flg) n=m+2;
44: PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL);
45: PetscPrintf(PETSC_COMM_WORLD,"\nRectangular bidiagonal matrix, m=%" PetscInt_FMT " n=%" PetscInt_FMT "\n\n",m,n);
47: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
48: Generate the matrix
49: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
51: MatCreate(PETSC_COMM_WORLD,&A);
52: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,n);
53: MatSetFromOptions(A);
54: MatSetUp(A);
55: MatGetOwnershipRange(A,&Istart,&Iend);
56: for (i=Istart;i<Iend;i++) {
57: col[0]=i; col[1]=i+1;
58: if (i<n-1) MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
59: else if (i==n-1) MatSetValue(A,i,col[0],value[0],INSERT_VALUES);
60: }
61: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
62: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
64: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
65: Compute singular values
66: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
68: SVDCreate(PETSC_COMM_WORLD,&svd);
69: SVDSetOperators(svd,A,NULL);
70: SVDSetType(svd,SVDPRIMME);
71: SVDPRIMMESetBlockSize(svd,4);
72: SVDPRIMMESetMethod(svd,SVD_PRIMME_HYBRID);
73: SVDSetFromOptions(svd);
75: SVDSolve(svd);
76: SVDPRIMMEGetBlockSize(svd,&bs);
77: SVDPRIMMEGetMethod(svd,&meth);
78: PetscPrintf(PETSC_COMM_WORLD," PRIMME: using block size %" PetscInt_FMT ", method %s\n",bs,SVDPRIMMEMethods[meth]);
80: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
81: Display solution and clean up
82: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
83: SVDErrorView(svd,SVD_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
84: SVDDestroy(&svd);
85: MatDestroy(&A);
86: SlepcFinalize();
87: return 0;
88: }
90: /*TEST
92: build:
93: requires: primme
95: test:
96: suffix: 1
97: nsize: 2
98: args: -svd_nsv 4
99: requires: primme
100: filter: sed -e "s/2.99255/2.99254/"
102: TEST*/