Line data Source code
1 : /*
2 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3 : SLEPc - Scalable Library for Eigenvalue Problem Computations
4 : Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
5 :
6 : This file is part of SLEPc.
7 : SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9 : */
10 :
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";
15 :
16 : #include <slepcsvd.h>
17 :
18 : /*
19 : This example computes the singular values of a rectangular bidiagonal matrix
20 :
21 : | 1 2 |
22 : | 1 2 |
23 : | 1 2 |
24 : A = | . . |
25 : | . . |
26 : | 1 2 |
27 : | 1 2 |
28 : */
29 :
30 2 : int main(int argc,char **argv)
31 : {
32 2 : Mat A;
33 2 : SVD svd;
34 2 : PetscInt m=20,n,Istart,Iend,i,k=6,col[2],bs;
35 2 : PetscScalar value[] = { 1, 2 };
36 2 : PetscBool flg;
37 2 : SVDPRIMMEMethod meth;
38 :
39 2 : PetscFunctionBeginUser;
40 2 : PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
41 2 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
42 2 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,&flg));
43 2 : if (!flg) n=m+2;
44 2 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL));
45 2 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nRectangular bidiagonal matrix, m=%" PetscInt_FMT " n=%" PetscInt_FMT "\n\n",m,n));
46 :
47 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
48 : Generate the matrix
49 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
50 :
51 2 : PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
52 2 : PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,n));
53 2 : PetscCall(MatSetFromOptions(A));
54 2 : PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
55 22 : for (i=Istart;i<Iend;i++) {
56 20 : col[0]=i; col[1]=i+1;
57 20 : if (i<n-1) PetscCall(MatSetValues(A,1,&i,2,col,value,INSERT_VALUES));
58 20 : else if (i==n-1) PetscCall(MatSetValue(A,i,col[0],value[0],INSERT_VALUES));
59 : }
60 2 : PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
61 2 : PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
62 :
63 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
64 : Compute singular values
65 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
66 :
67 2 : PetscCall(SVDCreate(PETSC_COMM_WORLD,&svd));
68 2 : PetscCall(SVDSetOperators(svd,A,NULL));
69 2 : PetscCall(SVDSetType(svd,SVDPRIMME));
70 2 : PetscCall(SVDPRIMMESetBlockSize(svd,4));
71 2 : PetscCall(SVDPRIMMESetMethod(svd,SVD_PRIMME_HYBRID));
72 2 : PetscCall(SVDSetFromOptions(svd));
73 :
74 2 : PetscCall(SVDSolve(svd));
75 2 : PetscCall(SVDPRIMMEGetBlockSize(svd,&bs));
76 2 : PetscCall(SVDPRIMMEGetMethod(svd,&meth));
77 2 : PetscCall(PetscPrintf(PETSC_COMM_WORLD," PRIMME: using block size %" PetscInt_FMT ", method %s\n",bs,SVDPRIMMEMethods[meth]));
78 :
79 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
80 : Display solution and clean up
81 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
82 2 : PetscCall(SVDErrorView(svd,SVD_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD));
83 2 : PetscCall(SVDDestroy(&svd));
84 2 : PetscCall(MatDestroy(&A));
85 2 : PetscCall(SlepcFinalize());
86 : return 0;
87 : }
88 :
89 : /*TEST
90 :
91 : build:
92 : requires: primme
93 :
94 : test:
95 : suffix: 1
96 : nsize: 2
97 : args: -svd_nsv 4
98 : requires: primme
99 : filter: sed -e "s/2.99255/2.99254/" | sed -e "s/2.97024/2.97023/"
100 :
101 : TEST*/
|