| Line | Branch | Exec | Source |
|---|---|---|---|
| 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 DSSVD.\n\n"; | ||
| 12 | |||
| 13 | #include <slepcds.h> | ||
| 14 | |||
| 15 | 40 | int main(int argc,char **argv) | |
| 16 | { | ||
| 17 | 40 | DS ds; | |
| 18 | 40 | SlepcSC sc; | |
| 19 | 40 | PetscReal sigma,rnorm,aux; | |
| 20 | 40 | PetscScalar *A,*U,*w,d; | |
| 21 | 40 | PetscInt i,j,k,n=15,m=10,m1,ld; | |
| 22 | 40 | PetscViewer viewer; | |
| 23 | 40 | PetscBool verbose,extrarow; | |
| 24 | |||
| 25 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
40 | PetscFunctionBeginUser; |
| 26 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
40 | PetscCall(SlepcInitialize(&argc,&argv,NULL,help)); |
| 27 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
40 | PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); |
| 28 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
40 | PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); |
| 29 | 40 | k = PetscMin(n,m); | |
| 30 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
40 | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type SVD - dimension %" PetscInt_FMT "x%" PetscInt_FMT ".\n",n,m)); |
| 31 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
40 | PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose)); |
| 32 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
40 | PetscCall(PetscOptionsHasName(NULL,NULL,"-extrarow",&extrarow)); |
| 33 | |||
| 34 | /* Create DS object */ | ||
| 35 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
40 | PetscCall(DSCreate(PETSC_COMM_WORLD,&ds)); |
| 36 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
40 | PetscCall(DSSetType(ds,DSSVD)); |
| 37 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
40 | PetscCall(DSSetFromOptions(ds)); |
| 38 | 40 | ld = PetscMax(n,m)+2; /* test leading dimension larger than n */ | |
| 39 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
40 | PetscCall(DSAllocate(ds,ld)); |
| 40 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
40 | PetscCall(DSSetDimensions(ds,n,0,0)); |
| 41 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
40 | PetscCall(DSSVDSetDimensions(ds,m)); |
| 42 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
40 | PetscCall(DSSVDGetDimensions(ds,&m1)); |
| 43 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
40 | PetscCheck(m1==m,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Inconsistent dimension value"); |
| 44 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
40 | PetscCall(DSSetExtraRow(ds,extrarow)); |
| 45 | |||
| 46 | /* Set up viewer */ | ||
| 47 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
40 | PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer)); |
| 48 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
40 | PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL)); |
| 49 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
40 | PetscCall(DSView(ds,viewer)); |
| 50 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
40 | PetscCall(PetscViewerPopFormat(viewer)); |
| 51 | |||
| 52 | /* Fill with a rectangular Toeplitz matrix */ | ||
| 53 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
40 | PetscCall(DSGetArray(ds,DS_MAT_A,&A)); |
| 54 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
440 | for (i=0;i<k;i++) A[i+i*ld]=1.0; |
| 55 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
120 | for (j=1;j<3;j++) { |
| 56 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
760 | for (i=0;i<m-j;i++) { |
| 57 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
680 | if ((i+j)<m) A[i+(i+j)*ld]=(PetscScalar)(j+1); |
| 58 | } | ||
| 59 | } | ||
| 60 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
280 | for (j=1;j<n/2;j++) { |
| 61 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
3000 | for (i=0;i<n-j;i++) { |
| 62 |
3/4✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
|
2760 | if ((i+j)<n && i<m) A[(i+j)+i*ld]=-1.0; |
| 63 | } | ||
| 64 | } | ||
| 65 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
40 | if (extrarow) { A[n-2+m*ld]=1.0; A[n-1+m*ld]=1.0; } /* really an extra column */ |
| 66 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
40 | PetscCall(DSRestoreArray(ds,DS_MAT_A,&A)); |
| 67 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
40 | PetscCall(DSSetState(ds,DS_STATE_RAW)); |
| 68 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
40 | if (verbose) { |
| 69 | ✗ | PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB)); | |
| 70 | ✗ | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n")); | |
| 71 | } | ||
| 72 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
40 | PetscCall(DSView(ds,viewer)); |
| 73 | |||
| 74 | /* Solve */ | ||
| 75 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
40 | PetscCall(PetscMalloc1(k,&w)); |
| 76 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
40 | PetscCall(DSGetSlepcSC(ds,&sc)); |
| 77 | 40 | sc->comparison = SlepcCompareLargestReal; | |
| 78 | 40 | sc->comparisonctx = NULL; | |
| 79 | 40 | sc->map = NULL; | |
| 80 | 40 | sc->mapobj = NULL; | |
| 81 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
40 | PetscCall(DSSolve(ds,w,NULL)); |
| 82 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
40 | PetscCall(DSSort(ds,w,NULL,NULL,NULL,NULL)); |
| 83 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
40 | if (extrarow) PetscCall(DSUpdateExtraRow(ds)); |
| 84 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
40 | if (verbose) { |
| 85 | ✗ | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n")); | |
| 86 | ✗ | PetscCall(DSView(ds,viewer)); | |
| 87 | } | ||
| 88 | |||
| 89 | /* Print singular values */ | ||
| 90 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
40 | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Computed singular values =\n")); |
| 91 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
440 | for (i=0;i<k;i++) { |
| 92 | 400 | sigma = PetscRealPart(w[i]); | |
| 93 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
400 | PetscCall(PetscViewerASCIIPrintf(viewer," %.5f\n",(double)sigma)); |
| 94 | } | ||
| 95 | |||
| 96 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
40 | if (extrarow) { |
| 97 | /* Check that extra column is correct */ | ||
| 98 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(DSGetArray(ds,DS_MAT_A,&A)); |
| 99 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(DSGetArray(ds,DS_MAT_U,&U)); |
| 100 | d = 0.0; | ||
| 101 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
320 | for (i=0;i<n;i++) d += A[i+m*ld]-U[n-2+i*ld]-U[n-1+i*ld]; |
| 102 |
1/8✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
20 | if (PetscAbsScalar(d)>10*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Warning: there is a mismatch in the extra row of %g\n",(double)PetscAbsScalar(d))); |
| 103 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(DSRestoreArray(ds,DS_MAT_A,&A)); |
| 104 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(DSRestoreArray(ds,DS_MAT_U,&U)); |
| 105 | } | ||
| 106 | |||
| 107 | /* Singular vectors */ | ||
| 108 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
40 | PetscCall(DSVectors(ds,DS_MAT_U,NULL,NULL)); /* all singular vectors */ |
| 109 | 40 | j = 0; | |
| 110 | 40 | rnorm = 0.0; | |
| 111 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
40 | PetscCall(DSGetArray(ds,DS_MAT_U,&U)); |
| 112 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
640 | for (i=0;i<n;i++) { |
| 113 | 600 | aux = PetscAbsScalar(U[i+j*ld]); | |
| 114 | 600 | rnorm += aux*aux; | |
| 115 | } | ||
| 116 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
40 | PetscCall(DSRestoreArray(ds,DS_MAT_U,&U)); |
| 117 | 40 | rnorm = PetscSqrtReal(rnorm); | |
| 118 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
40 | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Norm of 1st U vector = %.3f\n",(double)rnorm)); |
| 119 | |||
| 120 |
5/8✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
40 | PetscCall(PetscFree(w)); |
| 121 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
40 | PetscCall(DSDestroy(&ds)); |
| 122 |
3/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
40 | PetscCall(SlepcFinalize()); |
| 123 | return 0; | ||
| 124 | } | ||
| 125 | |||
| 126 | /*TEST | ||
| 127 | |||
| 128 | test: | ||
| 129 | args: -ds_method {{0 1}} | ||
| 130 | suffix: 1 | ||
| 131 | filter: grep -v "solving the problem" | ||
| 132 | requires: !single | ||
| 133 | |||
| 134 | test: | ||
| 135 | suffix: 2 | ||
| 136 | args: -extrarow -ds_method {{0 1}} | ||
| 137 | filter: grep -v "solving the problem" | ||
| 138 | requires: !single | ||
| 139 | |||
| 140 | TEST*/ | ||
| 141 |