| 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 | #include <slepc/private/dsimpl.h> /*I "slepcds.h" I*/ | ||
| 12 | #include <slepcblaslapack.h> | ||
| 13 | |||
| 14 | typedef struct { | ||
| 15 | PetscInt m; /* number of columns */ | ||
| 16 | PetscInt t; /* number of rows of V after truncating */ | ||
| 17 | } DS_SVD; | ||
| 18 | |||
| 19 | 1031 | static PetscErrorCode DSAllocate_SVD(DS ds,PetscInt ld) | |
| 20 | { | ||
| 21 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
1031 | PetscFunctionBegin; |
| 22 |
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.
|
1031 | if (!ds->compact) PetscCall(DSAllocateMat_Private(ds,DS_MAT_A)); |
| 23 |
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.
|
1031 | PetscCall(DSAllocateMat_Private(ds,DS_MAT_U)); |
| 24 |
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.
|
1031 | PetscCall(DSAllocateMat_Private(ds,DS_MAT_V)); |
| 25 |
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.
|
1031 | PetscCall(DSAllocateMat_Private(ds,DS_MAT_T)); |
| 26 |
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.
|
1031 | PetscCall(PetscFree(ds->perm)); |
| 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.
|
1031 | PetscCall(PetscMalloc1(ld,&ds->perm)); |
| 28 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
193 | PetscFunctionReturn(PETSC_SUCCESS); |
| 29 | } | ||
| 30 | |||
| 31 | /* 0 l k m-1 | ||
| 32 | ----------------------------------------- | ||
| 33 | |* . . | | ||
| 34 | | * . . | | ||
| 35 | | * . . | | ||
| 36 | | * . . | | ||
| 37 | | o o | | ||
| 38 | | o o | | ||
| 39 | | o o | | ||
| 40 | | o o | | ||
| 41 | | o o | | ||
| 42 | | o o | | ||
| 43 | | o x | | ||
| 44 | | x x | | ||
| 45 | | x x | | ||
| 46 | | x x | | ||
| 47 | | x x | | ||
| 48 | | x x | | ||
| 49 | | x x | | ||
| 50 | | x x | | ||
| 51 | | x x| | ||
| 52 | n-1 | x| | ||
| 53 | ----------------------------------------- | ||
| 54 | */ | ||
| 55 | |||
| 56 | 20 | static PetscErrorCode DSSwitchFormat_SVD(DS ds) | |
| 57 | { | ||
| 58 | 20 | DS_SVD *ctx = (DS_SVD*)ds->data; | |
| 59 | 20 | PetscReal *T; | |
| 60 | 20 | PetscScalar *A; | |
| 61 | 20 | PetscInt i,m=ctx->m,k=ds->k,ld=ds->ld; | |
| 62 | |||
| 63 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
20 | PetscFunctionBegin; |
| 64 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
20 | PetscCheck(m,PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"You should set the number of columns with DSSVDSetDimensions()"); |
| 65 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
20 | PetscCheck(ds->compact,PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Must have compact storage"); |
| 66 | /* switch from compact (arrow) to dense storage */ | ||
| 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.
|
20 | PetscCall(DSAllocateMat_Private(ds,DS_MAT_A)); |
| 68 |
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(MatDenseGetArrayWrite(ds->omat[DS_MAT_A],&A)); |
| 69 |
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(DSGetArrayReal(ds,DS_MAT_T,&T)); |
| 70 |
4/6✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(PetscArrayzero(A,ld*ld)); |
| 71 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
120 | for (i=0;i<k;i++) { |
| 72 | 100 | A[i+i*ld] = T[i]; | |
| 73 | 100 | A[i+k*ld] = T[i+ld]; | |
| 74 | } | ||
| 75 | 20 | A[k+k*ld] = T[k]; | |
| 76 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
100 | for (i=k+1;i<m;i++) { |
| 77 | 80 | A[i+i*ld] = T[i]; | |
| 78 | 80 | A[i-1+i*ld] = T[i-1+ld]; | |
| 79 | } | ||
| 80 |
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(MatDenseRestoreArrayWrite(ds->omat[DS_MAT_A],&A)); |
| 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.
|
20 | PetscCall(DSRestoreArrayReal(ds,DS_MAT_T,&T)); |
| 82 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
4 | PetscFunctionReturn(PETSC_SUCCESS); |
| 83 | } | ||
| 84 | |||
| 85 | 210 | static PetscErrorCode DSView_SVD(DS ds,PetscViewer viewer) | |
| 86 | { | ||
| 87 | 210 | DS_SVD *ctx = (DS_SVD*)ds->data; | |
| 88 | 210 | PetscViewerFormat format; | |
| 89 | 210 | PetscInt i,j,r,c,m=ctx->m,rows,cols; | |
| 90 | 210 | PetscReal *T,value; | |
| 91 | 210 | const char *methodname[] = { | |
| 92 | "Implicit zero-shift QR for bidiagonals (_bdsqr)", | ||
| 93 | "Divide and Conquer (_bdsdc or _gesdd)" | ||
| 94 | }; | ||
| 95 | 210 | const int nmeth=PETSC_STATIC_ARRAY_LENGTH(methodname); | |
| 96 | |||
| 97 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
210 | PetscFunctionBegin; |
| 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.
|
210 | PetscCall(PetscViewerGetFormat(viewer,&format)); |
| 99 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
210 | if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { |
| 100 |
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.
|
110 | PetscCall(PetscViewerASCIIPrintf(viewer,"number of columns: %" PetscInt_FMT "\n",m)); |
| 101 |
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.
|
110 | if (ds->method<nmeth) PetscCall(PetscViewerASCIIPrintf(viewer,"solving the problem with: %s\n",methodname[ds->method])); |
| 102 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
110 | PetscFunctionReturn(PETSC_SUCCESS); |
| 103 | } | ||
| 104 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
100 | PetscCheck(m,PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"You should set the number of columns with DSSVDSetDimensions()"); |
| 105 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
100 | if (ds->compact) { |
| 106 |
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.
|
60 | PetscCall(DSGetArrayReal(ds,DS_MAT_T,&T)); |
| 107 |
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.
|
60 | PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE)); |
| 108 | 60 | rows = ds->n; | |
| 109 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
60 | cols = ds->extrarow? m+1: m; |
| 110 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
60 | if (format == PETSC_VIEWER_ASCII_MATLAB) { |
| 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.
|
60 | PetscCall(PetscViewerASCIIPrintf(viewer,"%% Size = %" PetscInt_FMT " %" PetscInt_FMT "\n",rows,cols)); |
| 112 |
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.
|
60 | PetscCall(PetscViewerASCIIPrintf(viewer,"zzz = zeros(%" PetscInt_FMT ",3);\n",2*ds->n)); |
| 113 |
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.
|
60 | PetscCall(PetscViewerASCIIPrintf(viewer,"zzz = [\n")); |
| 114 |
7/8✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 2 times.
|
660 | for (i=0;i<PetscMin(ds->n,m);i++) PetscCall(PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT " %" PetscInt_FMT " %18.16e\n",i+1,i+1,(double)T[i])); |
| 115 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
620 | for (i=0;i<cols-1;i++) { |
| 116 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
560 | r = PetscMax(i+2,ds->k+1); |
| 117 | 560 | c = i+1; | |
| 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.
|
560 | PetscCall(PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT " %" PetscInt_FMT " %18.16e\n",c,r,(double)T[i+ds->ld])); |
| 119 | } | ||
| 120 |
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.
|
60 | PetscCall(PetscViewerASCIIPrintf(viewer,"];\n%s = spconvert(zzz);\n",DSMatName[DS_MAT_T])); |
| 121 | } else { | ||
| 122 | ✗ | for (i=0;i<rows;i++) { | |
| 123 | ✗ | for (j=0;j<cols;j++) { | |
| 124 | ✗ | if (i==j) value = T[i]; | |
| 125 | ✗ | else if (i<ds->k && j==ds->k) value = T[PetscMin(i,j)+ds->ld]; | |
| 126 | ✗ | else if (i+1==j && i>=ds->k) value = T[i+ds->ld]; | |
| 127 | else value = 0.0; | ||
| 128 | ✗ | PetscCall(PetscViewerASCIIPrintf(viewer," %18.16e ",(double)value)); | |
| 129 | } | ||
| 130 | ✗ | PetscCall(PetscViewerASCIIPrintf(viewer,"\n")); | |
| 131 | } | ||
| 132 | } | ||
| 133 |
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.
|
60 | PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE)); |
| 134 |
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.
|
60 | PetscCall(PetscViewerFlush(viewer)); |
| 135 |
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.
|
60 | PetscCall(DSRestoreArrayReal(ds,DS_MAT_T,&T)); |
| 136 |
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 | } else PetscCall(DSViewMat(ds,viewer,DS_MAT_A)); |
| 137 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
100 | if (ds->state>DS_STATE_INTERMEDIATE) { |
| 138 | ✗ | PetscCall(DSViewMat(ds,viewer,DS_MAT_U)); | |
| 139 | ✗ | PetscCall(DSViewMat(ds,viewer,DS_MAT_V)); | |
| 140 | } | ||
| 141 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
20 | PetscFunctionReturn(PETSC_SUCCESS); |
| 142 | } | ||
| 143 | |||
| 144 | 40 | static PetscErrorCode DSVectors_SVD(DS ds,DSMatType mat,PetscInt *j,PetscReal *rnorm) | |
| 145 | { | ||
| 146 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
40 | PetscFunctionBegin; |
| 147 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
40 | switch (mat) { |
| 148 | 40 | case DS_MAT_U: | |
| 149 | case DS_MAT_V: | ||
| 150 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
40 | if (rnorm) *rnorm = 0.0; |
| 151 | 8 | break; | |
| 152 | ✗ | default: | |
| 153 | ✗ | SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter"); | |
| 154 | } | ||
| 155 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
8 | PetscFunctionReturn(PETSC_SUCCESS); |
| 156 | } | ||
| 157 | |||
| 158 | 22319 | static PetscErrorCode DSSort_SVD(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k) | |
| 159 | { | ||
| 160 | 22319 | DS_SVD *ctx = (DS_SVD*)ds->data; | |
| 161 | 22319 | PetscInt n,l,i,*perm,ld=ds->ld; | |
| 162 | 22319 | PetscScalar *A; | |
| 163 | 22319 | PetscReal *d; | |
| 164 | |||
| 165 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
22319 | PetscFunctionBegin; |
| 166 |
2/14✓ Branch 0 taken 8 times.
✓ Branch 1 taken 2 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.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
|
22319 | if (!ds->sc) PetscFunctionReturn(PETSC_SUCCESS); |
| 167 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
22319 | PetscCheck(ctx->m,PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"You should set the number of columns with DSSVDSetDimensions()"); |
| 168 | 22319 | l = ds->l; | |
| 169 | 22319 | n = PetscMin(ds->n,ctx->m); | |
| 170 |
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.
|
22319 | PetscCall(DSGetArrayReal(ds,DS_MAT_T,&d)); |
| 171 | 22319 | perm = ds->perm; | |
| 172 |
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.
|
22319 | if (!rr) PetscCall(DSSortEigenvaluesReal_Private(ds,d,perm)); |
| 173 | ✗ | else PetscCall(DSSortEigenvalues_Private(ds,rr,ri,perm,PETSC_FALSE)); | |
| 174 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
243520 | for (i=l;i<n;i++) wr[i] = d[perm[i]]; |
| 175 |
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.
|
22319 | PetscCall(DSPermuteBoth_Private(ds,l,n,ds->n,ctx->m,DS_MAT_U,DS_MAT_V,perm)); |
| 176 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
243520 | for (i=l;i<n;i++) d[i] = PetscRealPart(wr[i]); |
| 177 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
22319 | if (!ds->compact) { |
| 178 |
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.
|
3568 | PetscCall(MatDenseGetArray(ds->omat[DS_MAT_A],&A)); |
| 179 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
42773 | for (i=l;i<n;i++) A[i+i*ld] = wr[i]; |
| 180 |
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.
|
3568 | PetscCall(MatDenseRestoreArray(ds->omat[DS_MAT_A],&A)); |
| 181 | } | ||
| 182 |
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.
|
22319 | PetscCall(DSRestoreArrayReal(ds,DS_MAT_T,&d)); |
| 183 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
3722 | PetscFunctionReturn(PETSC_SUCCESS); |
| 184 | } | ||
| 185 | |||
| 186 | 18731 | static PetscErrorCode DSUpdateExtraRow_SVD(DS ds) | |
| 187 | { | ||
| 188 | 18731 | DS_SVD *ctx = (DS_SVD*)ds->data; | |
| 189 | 18731 | PetscInt i; | |
| 190 | 18731 | PetscBLASInt n=0,m=0,ld,incx=1; | |
| 191 | 18731 | PetscScalar *A,*x,*y,one=1.0,zero=0.0; | |
| 192 | 18731 | PetscReal *T,*e,beta; | |
| 193 | 18731 | const PetscScalar *U; | |
| 194 | |||
| 195 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
18731 | PetscFunctionBegin; |
| 196 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
18731 | PetscCheck(ctx->m,PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"You should set the number of columns with DSSVDSetDimensions()"); |
| 197 |
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.
|
18731 | PetscCall(PetscBLASIntCast(ds->n,&n)); |
| 198 |
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.
|
18731 | PetscCall(PetscBLASIntCast(ctx->m,&m)); |
| 199 |
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.
|
18731 | PetscCall(PetscBLASIntCast(ds->ld,&ld)); |
| 200 |
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.
|
18731 | PetscCall(MatDenseGetArrayRead(ds->omat[DS_MAT_U],&U)); |
| 201 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
18731 | if (ds->compact) { |
| 202 |
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.
|
18711 | PetscCall(DSGetArrayReal(ds,DS_MAT_T,&T)); |
| 203 | 18711 | e = T+ld; | |
| 204 | 18711 | beta = e[m-1]; /* in compact, we assume all entries are zero except the last one */ | |
| 205 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
219121 | for (i=0;i<n;i++) e[i] = PetscRealPart(beta*U[n-1+i*ld]); |
| 206 | 18711 | ds->k = m; | |
| 207 |
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.
|
18711 | PetscCall(DSRestoreArrayReal(ds,DS_MAT_T,&T)); |
| 208 | } else { | ||
| 209 |
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(MatDenseGetArray(ds->omat[DS_MAT_A],&A)); |
| 210 |
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(DSAllocateWork_Private(ds,2*ld,0,0)); |
| 211 | 20 | x = ds->work; | |
| 212 | 20 | y = ds->work+ld; | |
| 213 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
320 | for (i=0;i<n;i++) x[i] = PetscConj(A[i+m*ld]); |
| 214 |
10/20✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
|
20 | PetscCallBLAS("BLASgemv",BLASgemv_("C",&n,&n,&one,U,&ld,x,&incx,&zero,y,&incx)); |
| 215 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
320 | for (i=0;i<n;i++) A[i+m*ld] = PetscConj(y[i]); |
| 216 | 20 | ds->k = m; | |
| 217 |
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(MatDenseRestoreArray(ds->omat[DS_MAT_A],&A)); |
| 218 | } | ||
| 219 |
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.
|
18731 | PetscCall(MatDenseRestoreArrayRead(ds->omat[DS_MAT_U],&U)); |
| 220 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
3062 | PetscFunctionReturn(PETSC_SUCCESS); |
| 221 | } | ||
| 222 | |||
| 223 | 12032 | static PetscErrorCode DSTruncate_SVD(DS ds,PetscInt n,PetscBool trim) | |
| 224 | { | ||
| 225 | 12032 | PetscInt i,ld=ds->ld,l=ds->l; | |
| 226 | 12032 | PetscScalar *A; | |
| 227 | 12032 | DS_SVD *ctx = (DS_SVD*)ds->data; | |
| 228 | |||
| 229 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
12032 | PetscFunctionBegin; |
| 230 |
1/10✗ 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.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
|
12032 | if (!ds->compact && ds->extrarow) PetscCall(MatDenseGetArray(ds->omat[DS_MAT_A],&A)); |
| 231 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
12032 | if (trim) { |
| 232 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
570 | if (!ds->compact && ds->extrarow) { /* clean extra column */ |
| 233 | ✗ | for (i=l;i<ds->n;i++) A[i+ctx->m*ld] = 0.0; | |
| 234 | } | ||
| 235 | 570 | ds->l = 0; | |
| 236 | 570 | ds->k = 0; | |
| 237 | 570 | ds->n = n; | |
| 238 | 570 | ctx->m = n; | |
| 239 | 570 | ds->t = ds->n; /* truncated length equal to the new dimension */ | |
| 240 | 570 | ctx->t = ctx->m; /* must also keep the previous dimension of V */ | |
| 241 | } else { | ||
| 242 |
1/6✗ 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.
|
11462 | if (!ds->compact && ds->extrarow && ds->k==ds->n) { |
| 243 | /* copy entries of extra column to the new position, then clean last row */ | ||
| 244 | ✗ | for (i=l;i<n;i++) A[i+n*ld] = A[i+ctx->m*ld]; | |
| 245 | ✗ | for (i=l;i<ds->n;i++) A[i+ctx->m*ld] = 0.0; | |
| 246 | } | ||
| 247 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
11462 | ds->k = (ds->extrarow)? n: 0; |
| 248 | 11462 | ds->t = ds->n; /* truncated length equal to previous dimension */ | |
| 249 | 11462 | ctx->t = ctx->m; /* must also keep the previous dimension of V */ | |
| 250 | 11462 | ds->n = n; | |
| 251 | 11462 | ctx->m = n; | |
| 252 | } | ||
| 253 |
1/10✗ 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.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
|
12032 | if (!ds->compact && ds->extrarow) PetscCall(MatDenseRestoreArray(ds->omat[DS_MAT_A],&A)); |
| 254 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
1776 | PetscFunctionReturn(PETSC_SUCCESS); |
| 255 | } | ||
| 256 | |||
| 257 | /* | ||
| 258 | DSArrowBidiag reduces a real square arrowhead matrix of the form | ||
| 259 | |||
| 260 | [ d 0 0 0 e ] | ||
| 261 | [ 0 d 0 0 e ] | ||
| 262 | A = [ 0 0 d 0 e ] | ||
| 263 | [ 0 0 0 d e ] | ||
| 264 | [ 0 0 0 0 d ] | ||
| 265 | |||
| 266 | to upper bidiagonal form | ||
| 267 | |||
| 268 | [ d e 0 0 0 ] | ||
| 269 | [ 0 d e 0 0 ] | ||
| 270 | B = Q'*A*P = [ 0 0 d e 0 ], | ||
| 271 | [ 0 0 0 d e ] | ||
| 272 | [ 0 0 0 0 d ] | ||
| 273 | |||
| 274 | where P,Q are orthogonal matrices. Uses plane rotations with a bulge chasing scheme. | ||
| 275 | On input, P and Q must be initialized to the identity matrix. | ||
| 276 | */ | ||
| 277 | 11482 | static PetscErrorCode DSArrowBidiag(PetscBLASInt n,PetscReal *d,PetscReal *e,PetscScalar *Q,PetscBLASInt ldq,PetscScalar *P,PetscBLASInt ldp) | |
| 278 | { | ||
| 279 | 11482 | PetscBLASInt i,j,j2,one=1; | |
| 280 | 11482 | PetscReal c,s,ct,st,off,temp0,temp1,temp2; | |
| 281 | |||
| 282 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
11482 | PetscFunctionBegin; |
| 283 |
2/14✓ Branch 0 taken 8 times.
✓ Branch 1 taken 2 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.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
|
11482 | if (n<=2) PetscFunctionReturn(PETSC_SUCCESS); |
| 284 | |||
| 285 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
54407 | for (j=0;j<n-2;j++) { |
| 286 | |||
| 287 | /* Eliminate entry e(j) by a rotation in the planes (j,j+1) */ | ||
| 288 | 42925 | temp0 = e[j+1]; | |
| 289 |
10/20✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
|
42925 | PetscCallBLAS("LAPACKlartg",LAPACKREALlartg_(&temp0,&e[j],&c,&s,&e[j+1])); |
| 290 | 42925 | s = -s; | |
| 291 | |||
| 292 | /* Apply rotation to Q */ | ||
| 293 | 42925 | j2 = j+2; | |
| 294 |
10/20✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
|
42925 | PetscCallBLAS("BLASrot",BLASMIXEDrot_(&j2,Q+j*ldq,&one,Q+(j+1)*ldq,&one,&c,&s)); |
| 295 | |||
| 296 | /* Apply rotation to diagonal elements, eliminate newly introduced entry A(j+1,j) */ | ||
| 297 | 42925 | temp0 = d[j+1]; | |
| 298 | 42925 | temp1 = c*temp0; | |
| 299 | 42925 | temp2 = -s*d[j]; | |
| 300 |
10/20✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
|
42925 | PetscCallBLAS("LAPACKlartg",LAPACKREALlartg_(&temp1,&temp2,&ct,&st,&d[j+1])); |
| 301 | 42925 | st = -st; | |
| 302 | 42925 | e[j] = -c*st*d[j] + s*ct*temp0; | |
| 303 | 42925 | d[j] = c*ct*d[j] + s*st*temp0; | |
| 304 | |||
| 305 | /* Apply rotation to P */ | ||
| 306 |
10/20✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
|
42925 | PetscCallBLAS("BLASrot",BLASMIXEDrot_(&j2,P+j*ldp,&one,P+(j+1)*ldp,&one,&ct,&st)); |
| 307 | |||
| 308 | /* Chase newly introduced off-diagonal entry to the top left corner */ | ||
| 309 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
111002 | for (i=j-1;i>=0;i--) { |
| 310 | |||
| 311 | /* Upper bulge */ | ||
| 312 | 68077 | off = -st*e[i]; | |
| 313 | 68077 | e[i] = ct*e[i]; | |
| 314 | 68077 | temp0 = e[i+1]; | |
| 315 |
10/20✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
|
68077 | PetscCallBLAS("LAPACKlartg",LAPACKREALlartg_(&temp0,&off,&c,&s,&e[i+1])); |
| 316 | 68077 | s = -s; | |
| 317 |
10/20✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
|
68077 | PetscCallBLAS("BLASrot",BLASMIXEDrot_(&j2,Q+i*ldq,&one,Q+(i+1)*ldq,&one,&c,&s)); |
| 318 | |||
| 319 | /* Lower bulge */ | ||
| 320 | 68077 | temp0 = d[i+1]; | |
| 321 | 68077 | temp1 = -s*e[i] + c*temp0; | |
| 322 | 68077 | temp2 = c*e[i] + s*temp0; | |
| 323 | 68077 | off = -s*d[i]; | |
| 324 |
10/20✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
|
68077 | PetscCallBLAS("LAPACKlartg",LAPACKREALlartg_(&temp1,&off,&ct,&st,&d[i+1])); |
| 325 | 68077 | st = -st; | |
| 326 | 68077 | e[i] = -c*st*d[i] + ct*temp2; | |
| 327 | 68077 | d[i] = c*ct*d[i] + st*temp2; | |
| 328 |
10/20✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
|
68077 | PetscCallBLAS("BLASrot",BLASMIXEDrot_(&j2,P+i*ldp,&one,P+(i+1)*ldp,&one,&ct,&st)); |
| 329 | } | ||
| 330 | } | ||
| 331 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
1678 | PetscFunctionReturn(PETSC_SUCCESS); |
| 332 | } | ||
| 333 | |||
| 334 | /* | ||
| 335 | Reduce to bidiagonal form by means of DSArrowBidiag. | ||
| 336 | */ | ||
| 337 | 23163 | static PetscErrorCode DSIntermediate_SVD(DS ds) | |
| 338 | { | ||
| 339 | 23163 | DS_SVD *ctx = (DS_SVD*)ds->data; | |
| 340 | 23163 | PetscInt i,j; | |
| 341 | 23163 | PetscBLASInt n1 = 0,n2,m2,lwork,info,l = 0,n = 0,m = 0,nm,ld,off; | |
| 342 | 23163 | PetscScalar *A,*U,*V,*W,*work,*tauq,*taup; | |
| 343 | 23163 | PetscReal *d,*e; | |
| 344 | |||
| 345 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
23163 | PetscFunctionBegin; |
| 346 |
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.
|
23163 | PetscCall(PetscBLASIntCast(ds->n,&n)); |
| 347 |
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.
|
23163 | PetscCall(PetscBLASIntCast(ctx->m,&m)); |
| 348 |
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.
|
23163 | PetscCall(PetscBLASIntCast(ds->l,&l)); |
| 349 |
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.
|
23163 | PetscCall(PetscBLASIntCast(ds->ld,&ld)); |
| 350 |
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.
|
23163 | PetscCall(PetscBLASIntCast(PetscMax(0,ds->k-l+1),&n1)); /* size of leading block, excl. locked */ |
| 351 | 23163 | n2 = n-l; /* n2 = n1 + size of trailing block */ | |
| 352 | 23163 | m2 = m-l; | |
| 353 | 23163 | off = l+l*ld; | |
| 354 | 23163 | nm = PetscMin(n,m); | |
| 355 |
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.
|
23163 | PetscCall(DSGetArrayReal(ds,DS_MAT_T,&d)); |
| 356 | 23163 | e = d+ld; | |
| 357 |
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.
|
23163 | PetscCall(MatDenseGetArray(ds->omat[DS_MAT_U],&U)); |
| 358 |
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.
|
23163 | PetscCall(MatDenseGetArray(ds->omat[DS_MAT_V],&V)); |
| 359 |
4/6✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
23163 | PetscCall(PetscArrayzero(U,ld*ld)); |
| 360 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
279024 | for (i=0;i<n;i++) U[i+i*ld] = 1.0; |
| 361 |
4/6✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
23163 | PetscCall(PetscArrayzero(V,ld*ld)); |
| 362 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
278744 | for (i=0;i<m;i++) V[i+i*ld] = 1.0; |
| 363 | |||
| 364 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
23163 | if (ds->compact) { |
| 365 | |||
| 366 |
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.
|
18721 | if (ds->state<DS_STATE_INTERMEDIATE) PetscCall(DSArrowBidiag(n1,d+l,e+l,U+off,ld,V+off,ld)); |
| 367 | |||
| 368 | } else { | ||
| 369 | |||
| 370 |
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.
|
4442 | PetscCall(MatDenseGetArray(ds->omat[DS_MAT_A],&A)); |
| 371 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
5391 | for (i=0;i<l;i++) { d[i] = PetscRealPart(A[i+i*ld]); e[i] = 0.0; } |
| 372 | |||
| 373 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
4442 | if (ds->state<DS_STATE_INTERMEDIATE) { |
| 374 | 4442 | lwork = (m+n)*16; | |
| 375 |
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.
|
4442 | PetscCall(DSAllocateWork_Private(ds,2*nm+ld*ld+lwork,0,0)); |
| 376 | 4442 | tauq = ds->work; | |
| 377 | 4442 | taup = ds->work+nm; | |
| 378 | 4442 | W = ds->work+2*nm; | |
| 379 | 4442 | work = ds->work+2*nm+ld*ld; | |
| 380 |
7/8✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 2 times.
|
59513 | for (j=0;j<m;j++) PetscCall(PetscArraycpy(W+j*ld,A+j*ld,n)); |
| 381 |
10/20✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
|
4442 | PetscCallBLAS("LAPACKgebrd",LAPACKgebrd_(&n2,&m2,W+off,&ld,d+l,e+l,tauq,taup,work,&lwork,&info)); |
| 382 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
4442 | SlepcCheckLapackInfo("gebrd",info); |
| 383 |
10/20✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
|
4442 | PetscCallBLAS("LAPACKormbr",LAPACKormbr_("Q","L","N",&n2,&n2,&m2,W+off,&ld,tauq,U+off,&ld,work,&lwork,&info)); |
| 384 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
4442 | SlepcCheckLapackInfo("ormbr",info); |
| 385 |
10/20✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
|
4442 | PetscCallBLAS("LAPACKormbr",LAPACKormbr_("P","R","N",&m2,&m2,&n2,W+off,&ld,taup,V+off,&ld,work,&lwork,&info)); |
| 386 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
4442 | SlepcCheckLapackInfo("ormbr",info); |
| 387 | } else { | ||
| 388 | /* copy bidiagonal to d,e */ | ||
| 389 | ✗ | for (i=l;i<nm;i++) d[i] = PetscRealPart(A[i+i*ld]); | |
| 390 | ✗ | for (i=l;i<nm-1;i++) e[i] = PetscRealPart(A[i+(i+1)*ld]); | |
| 391 | } | ||
| 392 |
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.
|
4442 | PetscCall(MatDenseRestoreArray(ds->omat[DS_MAT_A],&A)); |
| 393 | } | ||
| 394 |
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.
|
23163 | PetscCall(MatDenseRestoreArray(ds->omat[DS_MAT_U],&U)); |
| 395 |
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.
|
23163 | PetscCall(MatDenseRestoreArray(ds->omat[DS_MAT_V],&V)); |
| 396 |
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.
|
23163 | PetscCall(DSRestoreArrayReal(ds,DS_MAT_T,&d)); |
| 397 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
3836 | PetscFunctionReturn(PETSC_SUCCESS); |
| 398 | } | ||
| 399 | |||
| 400 | 23163 | static PetscErrorCode DSSolve_SVD_QR(DS ds,PetscScalar *wr,PetscScalar *wi) | |
| 401 | { | ||
| 402 | 23163 | DS_SVD *ctx = (DS_SVD*)ds->data; | |
| 403 | 23163 | PetscInt i,j; | |
| 404 | 23163 | PetscBLASInt n1,m1,info,l = 0,n = 0,m = 0,nm,ld,off,zero=0; | |
| 405 | 23163 | PetscScalar *A,*U,*V,*Vt; | |
| 406 | 23163 | PetscReal *d,*e; | |
| 407 | |||
| 408 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
23163 | PetscFunctionBegin; |
| 409 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
23163 | PetscCheck(ctx->m,PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"You should set the number of columns with DSSVDSetDimensions()"); |
| 410 |
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.
|
23163 | PetscCall(PetscBLASIntCast(ds->n,&n)); |
| 411 |
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.
|
23163 | PetscCall(PetscBLASIntCast(ctx->m,&m)); |
| 412 |
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.
|
23163 | PetscCall(PetscBLASIntCast(ds->l,&l)); |
| 413 |
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.
|
23163 | PetscCall(PetscBLASIntCast(ds->ld,&ld)); |
| 414 | 23163 | n1 = n-l; /* n1 = size of leading block, excl. locked + size of trailing block */ | |
| 415 | 23163 | m1 = m-l; | |
| 416 | 23163 | nm = PetscMin(n1,m1); | |
| 417 | 23163 | off = l+l*ld; | |
| 418 |
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.
|
23163 | PetscCall(DSGetArrayReal(ds,DS_MAT_T,&d)); |
| 419 | 23163 | e = d+ld; | |
| 420 | |||
| 421 | /* Reduce to bidiagonal form */ | ||
| 422 |
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.
|
23163 | PetscCall(DSIntermediate_SVD(ds)); |
| 423 | |||
| 424 |
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.
|
23163 | PetscCall(MatDenseGetArray(ds->omat[DS_MAT_U],&U)); |
| 425 |
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.
|
23163 | PetscCall(MatDenseGetArray(ds->omat[DS_MAT_V],&V)); |
| 426 | |||
| 427 | /* solve bidiagonal SVD problem */ | ||
| 428 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
42886 | for (i=0;i<l;i++) wr[i] = d[i]; |
| 429 |
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.
|
23163 | PetscCall(DSAllocateWork_Private(ds,ld*ld,4*n1,0)); |
| 430 | 23163 | Vt = ds->work; | |
| 431 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
259021 | for (i=l;i<m;i++) { |
| 432 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2949816 | for (j=l;j<m;j++) { |
| 433 | 2713958 | Vt[i+j*ld] = PetscConj(V[j+i*ld]); /* Lapack expects transposed VT */ | |
| 434 | } | ||
| 435 | } | ||
| 436 |
13/22✓ Branch 0 taken 10 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 2 times.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
✓ Branch 19 taken 2 times.
✗ Branch 20 not taken.
✗ Branch 21 not taken.
|
23193 | PetscCallBLAS("LAPACKbdsqr",LAPACKbdsqr_(n>=m?"U":"L",&nm,&m1,&n1,&zero,d+l,e+l,Vt+off,&ld,U+off,&ld,NULL,&ld,ds->rwork,&info)); |
| 437 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
23163 | SlepcCheckLapackInfo("bdsqr",info); |
| 438 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
259021 | for (i=l;i<m;i++) { |
| 439 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2949816 | for (j=l;j<m;j++) { |
| 440 | 2713958 | V[i+j*ld] = PetscConj(Vt[j+i*ld]); /* transpose VT returned by Lapack */ | |
| 441 | } | ||
| 442 | } | ||
| 443 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
258961 | for (i=l;i<PetscMin(ds->n,ctx->m);i++) wr[i] = d[i]; |
| 444 | |||
| 445 | /* create diagonal matrix as a result */ | ||
| 446 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 5 times.
✓ Branch 3 taken 5 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
23163 | if (ds->compact) PetscCall(PetscArrayzero(e,n-1)); |
| 447 | else { | ||
| 448 |
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.
|
4442 | PetscCall(MatDenseGetArray(ds->omat[DS_MAT_A],&A)); |
| 449 |
7/8✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 2 times.
|
58564 | for (i=l;i<m;i++) PetscCall(PetscArrayzero(A+l+i*ld,n-l)); |
| 450 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
58504 | for (i=l;i<PetscMin(ds->n,ctx->m);i++) A[i+i*ld] = d[i]; |
| 451 |
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.
|
4442 | PetscCall(MatDenseRestoreArray(ds->omat[DS_MAT_A],&A)); |
| 452 | } | ||
| 453 |
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.
|
23163 | PetscCall(MatDenseRestoreArray(ds->omat[DS_MAT_U],&U)); |
| 454 |
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.
|
23163 | PetscCall(MatDenseRestoreArray(ds->omat[DS_MAT_V],&V)); |
| 455 |
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.
|
23163 | PetscCall(DSRestoreArrayReal(ds,DS_MAT_T,&d)); |
| 456 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
3836 | PetscFunctionReturn(PETSC_SUCCESS); |
| 457 | } | ||
| 458 | |||
| 459 | 50 | static PetscErrorCode DSSolve_SVD_DC(DS ds,PetscScalar *wr,PetscScalar *wi) | |
| 460 | { | ||
| 461 | 50 | DS_SVD *ctx = (DS_SVD*)ds->data; | |
| 462 | 50 | PetscInt i,j; | |
| 463 | 50 | PetscBLASInt n1,m1,info,l = 0,n = 0,m = 0,nm,ld,off,lwork; | |
| 464 | 50 | PetscScalar *A,*U,*V,*W,qwork; | |
| 465 | 50 | PetscReal *d,*e,*Ur,*Vr; | |
| 466 | |||
| 467 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
50 | PetscFunctionBegin; |
| 468 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
50 | PetscCheck(ctx->m,PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"You should set the number of columns with DSSVDSetDimensions()"); |
| 469 |
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.
|
50 | PetscCall(PetscBLASIntCast(ds->n,&n)); |
| 470 |
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.
|
50 | PetscCall(PetscBLASIntCast(ctx->m,&m)); |
| 471 |
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.
|
50 | PetscCall(PetscBLASIntCast(ds->l,&l)); |
| 472 |
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.
|
50 | PetscCall(PetscBLASIntCast(ds->ld,&ld)); |
| 473 | 50 | n1 = n-l; /* n1 = size of leading block, excl. locked + size of trailing block */ | |
| 474 | 50 | m1 = m-l; | |
| 475 | 50 | off = l+l*ld; | |
| 476 |
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.
|
50 | if (ds->compact) PetscCall(DSAllocateMat_Private(ds,DS_MAT_A)); |
| 477 |
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.
|
50 | PetscCall(MatDenseGetArray(ds->omat[DS_MAT_A],&A)); |
| 478 |
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.
|
50 | PetscCall(MatDenseGetArrayWrite(ds->omat[DS_MAT_U],&U)); |
| 479 |
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.
|
50 | PetscCall(MatDenseGetArrayWrite(ds->omat[DS_MAT_V],&V)); |
| 480 |
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.
|
50 | PetscCall(DSGetArrayReal(ds,DS_MAT_T,&d)); |
| 481 | 50 | e = d+ld; | |
| 482 |
4/6✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
50 | PetscCall(PetscArrayzero(U,ld*ld)); |
| 483 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
90 | for (i=0;i<l;i++) U[i+i*ld] = 1.0; |
| 484 |
4/6✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
50 | PetscCall(PetscArrayzero(V,ld*ld)); |
| 485 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
90 | for (i=0;i<l;i++) V[i+i*ld] = 1.0; |
| 486 | |||
| 487 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
50 | if (ds->state>DS_STATE_RAW) { |
| 488 | /* solve bidiagonal SVD problem */ | ||
| 489 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
10 | for (i=0;i<l;i++) wr[i] = d[i]; |
| 490 | #if defined(PETSC_USE_COMPLEX) | ||
| 491 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(DSAllocateWork_Private(ds,0,3*n1*n1+4*n1+2*ld*ld,8*n1)); |
| 492 | 5 | Ur = ds->rwork+3*n1*n1+4*n1; | |
| 493 | 5 | Vr = ds->rwork+3*n1*n1+4*n1+ld*ld; | |
| 494 | #else | ||
| 495 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(DSAllocateWork_Private(ds,0,3*n1*n1+4*n1+ld*ld,8*n1)); |
| 496 | 5 | Ur = U; | |
| 497 | 5 | Vr = ds->rwork+3*n1*n1+4*n1; | |
| 498 | #endif | ||
| 499 |
10/20✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
|
10 | PetscCallBLAS("LAPACKbdsdc",LAPACKbdsdc_("U","I",&n1,d+l,e+l,Ur+off,&ld,Vr+off,&ld,NULL,NULL,ds->rwork,ds->iwork,&info)); |
| 500 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
10 | SlepcCheckLapackInfo("bdsdc",info); |
| 501 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
110 | for (i=l;i<n;i++) { |
| 502 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1100 | for (j=l;j<n;j++) { |
| 503 | #if defined(PETSC_USE_COMPLEX) | ||
| 504 | 500 | U[i+j*ld] = Ur[i+j*ld]; | |
| 505 | #endif | ||
| 506 | 1000 | V[i+j*ld] = PetscConj(Vr[j+i*ld]); /* transpose VT returned by Lapack */ | |
| 507 | } | ||
| 508 | } | ||
| 509 | } else { | ||
| 510 | /* solve general rectangular SVD problem */ | ||
| 511 |
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(DSAllocateMat_Private(ds,DS_MAT_W)); |
| 512 |
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(MatDenseGetArrayWrite(ds->omat[DS_MAT_W],&W)); |
| 513 |
5/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
|
40 | if (ds->compact) PetscCall(DSSwitchFormat_SVD(ds)); |
| 514 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
80 | for (i=0;i<l;i++) wr[i] = d[i]; |
| 515 | 40 | nm = PetscMin(n,m); | |
| 516 |
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(DSAllocateWork_Private(ds,0,0,8*nm)); |
| 517 | 40 | lwork = -1; | |
| 518 | #if defined(PETSC_USE_COMPLEX) | ||
| 519 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
20 | PetscCall(DSAllocateWork_Private(ds,0,5*nm*nm+7*nm,0)); |
| 520 |
10/20✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 1 times.
✓ Branch 8 taken 1 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 1 times.
✓ Branch 12 taken 1 times.
✗ Branch 13 not taken.
✓ Branch 14 taken 1 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 1 times.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
|
20 | PetscCallBLAS("LAPACKgesdd",LAPACKgesdd_("A",&n1,&m1,A+off,&ld,d+l,U+off,&ld,W+off,&ld,&qwork,&lwork,ds->rwork,ds->iwork,&info)); |
| 521 | #else | ||
| 522 |
10/20✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 1 times.
✓ Branch 8 taken 1 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 1 times.
✓ Branch 12 taken 1 times.
✗ Branch 13 not taken.
✓ Branch 14 taken 1 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 1 times.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
|
20 | PetscCallBLAS("LAPACKgesdd",LAPACKgesdd_("A",&n1,&m1,A+off,&ld,d+l,U+off,&ld,W+off,&ld,&qwork,&lwork,ds->iwork,&info)); |
| 523 | #endif | ||
| 524 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
40 | SlepcCheckLapackInfo("gesdd",info); |
| 525 |
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(PetscBLASIntCast((PetscInt)PetscRealPart(qwork),&lwork)); |
| 526 |
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(DSAllocateWork_Private(ds,lwork,0,0)); |
| 527 | #if defined(PETSC_USE_COMPLEX) | ||
| 528 |
10/20✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 1 times.
✓ Branch 8 taken 1 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 1 times.
✓ Branch 12 taken 1 times.
✗ Branch 13 not taken.
✓ Branch 14 taken 1 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 1 times.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
|
20 | PetscCallBLAS("LAPACKgesdd",LAPACKgesdd_("A",&n1,&m1,A+off,&ld,d+l,U+off,&ld,W+off,&ld,ds->work,&lwork,ds->rwork,ds->iwork,&info)); |
| 529 | #else | ||
| 530 |
10/20✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 1 times.
✓ Branch 8 taken 1 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 1 times.
✓ Branch 12 taken 1 times.
✗ Branch 13 not taken.
✓ Branch 14 taken 1 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 1 times.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
|
20 | PetscCallBLAS("LAPACKgesdd",LAPACKgesdd_("A",&n1,&m1,A+off,&ld,d+l,U+off,&ld,W+off,&ld,ds->work,&lwork,ds->iwork,&info)); |
| 531 | #endif | ||
| 532 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
40 | SlepcCheckLapackInfo("gesdd",info); |
| 533 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
400 | for (i=l;i<m;i++) { |
| 534 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
3640 | for (j=l;j<m;j++) V[i+j*ld] = PetscConj(W[j+i*ld]); /* transpose VT returned by Lapack */ |
| 535 | } | ||
| 536 |
5/6✓ Branch 0 taken 9 times.
✓ Branch 1 taken 1 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 1 times.
✓ Branch 5 taken 1 times.
|
40 | PetscCall(MatDenseRestoreArrayWrite(ds->omat[DS_MAT_W],&W)); |
| 537 | } | ||
| 538 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
510 | for (i=l;i<PetscMin(ds->n,ctx->m);i++) wr[i] = d[i]; |
| 539 | |||
| 540 | /* create diagonal matrix as a result */ | ||
| 541 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 5 times.
✓ Branch 3 taken 5 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
50 | if (ds->compact) PetscCall(PetscArrayzero(e,n-1)); |
| 542 | else { | ||
| 543 |
7/8✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 2 times.
|
220 | for (i=l;i<m;i++) PetscCall(PetscArrayzero(A+l+i*ld,n-l)); |
| 544 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
320 | for (i=l;i<n;i++) A[i+i*ld] = d[i]; |
| 545 | } | ||
| 546 |
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.
|
50 | PetscCall(MatDenseRestoreArray(ds->omat[DS_MAT_A],&A)); |
| 547 |
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.
|
50 | PetscCall(MatDenseRestoreArrayWrite(ds->omat[DS_MAT_U],&U)); |
| 548 |
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.
|
50 | PetscCall(MatDenseRestoreArrayWrite(ds->omat[DS_MAT_V],&V)); |
| 549 |
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.
|
50 | PetscCall(DSRestoreArrayReal(ds,DS_MAT_T,&d)); |
| 550 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
10 | PetscFunctionReturn(PETSC_SUCCESS); |
| 551 | } | ||
| 552 | |||
| 553 | #if !defined(PETSC_HAVE_MPIUNI) | ||
| 554 | 60 | static PetscErrorCode DSSynchronize_SVD(DS ds,PetscScalar eigr[],PetscScalar eigi[]) | |
| 555 | { | ||
| 556 | 60 | PetscInt ld=ds->ld,l=ds->l,k=0,kr=0; | |
| 557 | 60 | PetscMPIInt n,rank,off=0,size,ldn,ld3; | |
| 558 | 60 | PetscScalar *A,*U,*V; | |
| 559 | 60 | PetscReal *T; | |
| 560 | |||
| 561 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
60 | PetscFunctionBegin; |
| 562 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
60 | if (ds->compact) kr = 3*ld; |
| 563 | ✗ | else k = (ds->n-l)*ld; | |
| 564 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
60 | if (ds->state>DS_STATE_RAW) k += 2*(ds->n-l)*ld; |
| 565 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
60 | if (eigr) k += ds->n-l; |
| 566 |
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.
|
60 | PetscCall(DSAllocateWork_Private(ds,k+kr,0,0)); |
| 567 |
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.
|
60 | PetscCall(PetscMPIIntCast(k*sizeof(PetscScalar)+kr*sizeof(PetscReal),&size)); |
| 568 |
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.
|
60 | PetscCall(PetscMPIIntCast(ds->n-l,&n)); |
| 569 |
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.
|
60 | PetscCall(PetscMPIIntCast(ld*(ds->n-l),&ldn)); |
| 570 |
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.
|
60 | PetscCall(PetscMPIIntCast(3*ld,&ld3)); |
| 571 |
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.
|
60 | if (ds->compact) PetscCall(DSGetArrayReal(ds,DS_MAT_T,&T)); |
| 572 | ✗ | else PetscCall(MatDenseGetArray(ds->omat[DS_MAT_A],&A)); | |
| 573 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
60 | if (ds->state>DS_STATE_RAW) { |
| 574 |
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.
|
60 | PetscCall(MatDenseGetArray(ds->omat[DS_MAT_U],&U)); |
| 575 |
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.
|
60 | PetscCall(MatDenseGetArray(ds->omat[DS_MAT_V],&V)); |
| 576 | } | ||
| 577 |
14/28✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✓ Branch 21 taken 2 times.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
✗ Branch 24 not taken.
✓ Branch 25 taken 2 times.
✗ Branch 26 not taken.
✗ Branch 27 not taken.
|
60 | PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)ds),&rank)); |
| 578 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
60 | if (!rank) { |
| 579 |
15/30✓ 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 taken 8 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✗ Branch 24 not taken.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✗ Branch 28 not taken.
✗ Branch 29 not taken.
|
30 | if (ds->compact) PetscCallMPI(MPI_Pack(T,ld3,MPIU_REAL,ds->work,size,&off,PetscObjectComm((PetscObject)ds))); |
| 580 | ✗ | else PetscCallMPI(MPI_Pack(A+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds))); | |
| 581 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
30 | if (ds->state>DS_STATE_RAW) { |
| 582 |
14/28✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✓ Branch 21 taken 2 times.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
✗ Branch 24 not taken.
✓ Branch 25 taken 2 times.
✗ Branch 26 not taken.
✗ Branch 27 not taken.
|
30 | PetscCallMPI(MPI_Pack(U+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds))); |
| 583 |
14/28✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✓ Branch 21 taken 2 times.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
✗ Branch 24 not taken.
✓ Branch 25 taken 2 times.
✗ Branch 26 not taken.
✗ Branch 27 not taken.
|
30 | PetscCallMPI(MPI_Pack(V+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds))); |
| 584 | } | ||
| 585 |
15/30✓ 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 taken 8 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✗ Branch 24 not taken.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✗ Branch 28 not taken.
✗ Branch 29 not taken.
|
30 | if (eigr) PetscCallMPI(MPI_Pack(eigr+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds))); |
| 586 | } | ||
| 587 |
15/30✓ 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 taken 8 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✗ Branch 24 not taken.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✗ Branch 28 not taken.
✗ Branch 29 not taken.
|
120 | PetscCallMPI(MPI_Bcast(ds->work,size,MPI_BYTE,0,PetscObjectComm((PetscObject)ds))); |
| 588 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
60 | if (rank) { |
| 589 |
15/30✓ 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 taken 8 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✗ Branch 24 not taken.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✗ Branch 28 not taken.
✗ Branch 29 not taken.
|
30 | if (ds->compact) PetscCallMPI(MPI_Unpack(ds->work,size,&off,T,ld3,MPIU_REAL,PetscObjectComm((PetscObject)ds))); |
| 590 | ✗ | else PetscCallMPI(MPI_Unpack(ds->work,size,&off,A+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds))); | |
| 591 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
30 | if (ds->state>DS_STATE_RAW) { |
| 592 |
14/28✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✓ Branch 21 taken 2 times.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
✗ Branch 24 not taken.
✓ Branch 25 taken 2 times.
✗ Branch 26 not taken.
✗ Branch 27 not taken.
|
30 | PetscCallMPI(MPI_Unpack(ds->work,size,&off,U+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds))); |
| 593 |
14/28✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✓ Branch 21 taken 2 times.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
✗ Branch 24 not taken.
✓ Branch 25 taken 2 times.
✗ Branch 26 not taken.
✗ Branch 27 not taken.
|
30 | PetscCallMPI(MPI_Unpack(ds->work,size,&off,V+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds))); |
| 594 | } | ||
| 595 |
15/30✓ 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 taken 8 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✗ Branch 24 not taken.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✗ Branch 28 not taken.
✗ Branch 29 not taken.
|
30 | if (eigr) PetscCallMPI(MPI_Unpack(ds->work,size,&off,eigr+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds))); |
| 596 | } | ||
| 597 |
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.
|
60 | if (ds->compact) PetscCall(DSRestoreArrayReal(ds,DS_MAT_T,&T)); |
| 598 | ✗ | else PetscCall(MatDenseRestoreArray(ds->omat[DS_MAT_A],&A)); | |
| 599 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
60 | if (ds->state>DS_STATE_RAW) { |
| 600 |
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.
|
60 | PetscCall(MatDenseRestoreArray(ds->omat[DS_MAT_U],&U)); |
| 601 |
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.
|
60 | PetscCall(MatDenseRestoreArray(ds->omat[DS_MAT_V],&V)); |
| 602 | } | ||
| 603 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
12 | PetscFunctionReturn(PETSC_SUCCESS); |
| 604 | } | ||
| 605 | #endif | ||
| 606 | |||
| 607 | 48794 | static PetscErrorCode DSMatGetSize_SVD(DS ds,DSMatType t,PetscInt *rows,PetscInt *cols) | |
| 608 | { | ||
| 609 | 48794 | DS_SVD *ctx = (DS_SVD*)ds->data; | |
| 610 | |||
| 611 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
48794 | PetscFunctionBegin; |
| 612 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
48794 | PetscCheck(ctx->m,PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"You should set the number of columns with DSSVDSetDimensions()"); |
| 613 |
3/5✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
✗ Branch 4 not taken.
|
48794 | switch (t) { |
| 614 | 4462 | case DS_MAT_A: | |
| 615 | 4462 | *rows = ds->n; | |
| 616 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
4462 | *cols = ds->extrarow? ctx->m+1: ctx->m; |
| 617 | 4462 | break; | |
| 618 | ✗ | case DS_MAT_T: | |
| 619 | ✗ | *rows = ds->n; | |
| 620 | ✗ | *cols = PetscDefined(USE_COMPLEX)? 2: 3; | |
| 621 | ✗ | break; | |
| 622 | 22269 | case DS_MAT_U: | |
| 623 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
22269 | *rows = ds->state==DS_STATE_TRUNCATED? ds->t: ds->n; |
| 624 | 22269 | *cols = ds->n; | |
| 625 | 22269 | break; | |
| 626 | 22063 | case DS_MAT_V: | |
| 627 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
22063 | *rows = ds->state==DS_STATE_TRUNCATED? ctx->t: ctx->m; |
| 628 | 22063 | *cols = ctx->m; | |
| 629 | 22063 | break; | |
| 630 | ✗ | default: | |
| 631 | ✗ | SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid t parameter"); | |
| 632 | } | ||
| 633 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
8136 | PetscFunctionReturn(PETSC_SUCCESS); |
| 634 | } | ||
| 635 | |||
| 636 | 23213 | static PetscErrorCode DSSVDSetDimensions_SVD(DS ds,PetscInt m) | |
| 637 | { | ||
| 638 | 23213 | DS_SVD *ctx = (DS_SVD*)ds->data; | |
| 639 | |||
| 640 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
23213 | PetscFunctionBegin; |
| 641 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
23213 | DSCheckAlloc(ds,1); |
| 642 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
23213 | if (m==PETSC_DECIDE || m==PETSC_DEFAULT) { |
| 643 | ✗ | ctx->m = ds->ld; | |
| 644 | } else { | ||
| 645 |
2/6✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 10 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
23213 | PetscCheck(m>0 && m<=ds->ld,PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of m. Must be between 1 and ld"); |
| 646 | 23213 | ctx->m = m; | |
| 647 | } | ||
| 648 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
3846 | PetscFunctionReturn(PETSC_SUCCESS); |
| 649 | } | ||
| 650 | |||
| 651 | /*@ | ||
| 652 | DSSVDSetDimensions - Sets the number of columns for a `DSSVD`. | ||
| 653 | |||
| 654 | Logically Collective | ||
| 655 | |||
| 656 | Input Parameters: | ||
| 657 | + ds - the direct solver context | ||
| 658 | - m - the number of columns | ||
| 659 | |||
| 660 | Notes: | ||
| 661 | This call is complementary to `DSSetDimensions()`, to provide a dimension | ||
| 662 | that is specific to this `DS` type. | ||
| 663 | |||
| 664 | Level: intermediate | ||
| 665 | |||
| 666 | .seealso: [](sec:ds), `DSSVD`, `DSSVDGetDimensions()`, `DSSetDimensions()` | ||
| 667 | @*/ | ||
| 668 | 23213 | PetscErrorCode DSSVDSetDimensions(DS ds,PetscInt m) | |
| 669 | { | ||
| 670 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
23213 | PetscFunctionBegin; |
| 671 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
23213 | PetscValidHeaderSpecific(ds,DS_CLASSID,1); |
| 672 |
27/62✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
✓ Branch 19 taken 2 times.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✓ Branch 24 taken 2 times.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✓ Branch 28 taken 2 times.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
✓ Branch 31 taken 2 times.
✓ Branch 32 taken 2 times.
✗ Branch 33 not taken.
✓ Branch 34 taken 2 times.
✗ Branch 35 not taken.
✗ Branch 36 not taken.
✓ Branch 37 taken 2 times.
✗ Branch 38 not taken.
✗ Branch 39 not taken.
✗ Branch 40 not taken.
✓ Branch 41 taken 2 times.
✗ Branch 42 not taken.
✗ Branch 43 not taken.
✗ Branch 44 not taken.
✓ Branch 45 taken 2 times.
✓ Branch 46 taken 2 times.
✗ Branch 47 not taken.
✗ Branch 48 not taken.
✓ Branch 49 taken 2 times.
✓ Branch 50 taken 2 times.
✗ Branch 51 not taken.
✓ Branch 52 taken 2 times.
✗ Branch 53 not taken.
✗ Branch 54 not taken.
✓ Branch 55 taken 2 times.
✗ Branch 56 not taken.
✗ Branch 57 not taken.
✗ Branch 58 not taken.
✓ Branch 59 taken 2 times.
✗ Branch 60 not taken.
✗ Branch 61 not taken.
|
23213 | PetscValidLogicalCollectiveInt(ds,m,2); |
| 673 |
8/14✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 10 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
|
23213 | PetscTryMethod(ds,"DSSVDSetDimensions_C",(DS,PetscInt),(ds,m)); |
| 674 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
23213 | PetscFunctionReturn(PETSC_SUCCESS); |
| 675 | } | ||
| 676 | |||
| 677 | 40 | static PetscErrorCode DSSVDGetDimensions_SVD(DS ds,PetscInt *m) | |
| 678 | { | ||
| 679 | 40 | DS_SVD *ctx = (DS_SVD*)ds->data; | |
| 680 | |||
| 681 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
40 | PetscFunctionBegin; |
| 682 | 40 | *m = ctx->m; | |
| 683 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
40 | PetscFunctionReturn(PETSC_SUCCESS); |
| 684 | } | ||
| 685 | |||
| 686 | /*@ | ||
| 687 | DSSVDGetDimensions - Returns the number of columns for a `DSSVD`. | ||
| 688 | |||
| 689 | Not Collective | ||
| 690 | |||
| 691 | Input Parameter: | ||
| 692 | . ds - the direct solver context | ||
| 693 | |||
| 694 | Output Parameter: | ||
| 695 | . m - the number of columns | ||
| 696 | |||
| 697 | Level: intermediate | ||
| 698 | |||
| 699 | .seealso: [](sec:ds), `DSSVD`, `DSSVDSetDimensions()` | ||
| 700 | @*/ | ||
| 701 | 40 | PetscErrorCode DSSVDGetDimensions(DS ds,PetscInt *m) | |
| 702 | { | ||
| 703 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
40 | PetscFunctionBegin; |
| 704 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
40 | PetscValidHeaderSpecific(ds,DS_CLASSID,1); |
| 705 |
2/8✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
40 | PetscAssertPointer(m,2); |
| 706 |
9/16✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 10 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
|
40 | PetscUseMethod(ds,"DSSVDGetDimensions_C",(DS,PetscInt*),(ds,m)); |
| 707 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
40 | PetscFunctionReturn(PETSC_SUCCESS); |
| 708 | } | ||
| 709 | |||
| 710 | 1001 | static PetscErrorCode DSDestroy_SVD(DS ds) | |
| 711 | { | ||
| 712 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
1001 | PetscFunctionBegin; |
| 713 |
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.
|
1001 | PetscCall(PetscFree(ds->data)); |
| 714 |
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.
|
1001 | PetscCall(PetscObjectComposeFunction((PetscObject)ds,"DSSVDSetDimensions_C",NULL)); |
| 715 |
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.
|
1001 | PetscCall(PetscObjectComposeFunction((PetscObject)ds,"DSSVDGetDimensions_C",NULL)); |
| 716 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
187 | PetscFunctionReturn(PETSC_SUCCESS); |
| 717 | } | ||
| 718 | |||
| 719 | 60 | static PetscErrorCode DSSetCompact_SVD(DS ds,PetscBool comp) | |
| 720 | { | ||
| 721 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
60 | PetscFunctionBegin; |
| 722 |
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.
|
60 | if (!comp) PetscCall(DSAllocateMat_Private(ds,DS_MAT_A)); |
| 723 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
12 | PetscFunctionReturn(PETSC_SUCCESS); |
| 724 | } | ||
| 725 | |||
| 726 | 31 | static PetscErrorCode DSReallocate_SVD(DS ds,PetscInt ld) | |
| 727 | { | ||
| 728 | 31 | PetscInt i,*perm=ds->perm; | |
| 729 | |||
| 730 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
31 | PetscFunctionBegin; |
| 731 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
713 | for (i=0;i<DS_NUM_MAT;i++) { |
| 732 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
682 | if (!ds->compact && i==DS_MAT_A) continue; |
| 733 |
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.
|
682 | if (i!=DS_MAT_U && i!=DS_MAT_V && i!=DS_MAT_T) PetscCall(MatDestroy(&ds->omat[i])); |
| 734 | } | ||
| 735 | |||
| 736 |
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.
|
31 | if (!ds->compact) PetscCall(DSReallocateMat_Private(ds,DS_MAT_A,ld)); |
| 737 |
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.
|
31 | PetscCall(DSReallocateMat_Private(ds,DS_MAT_U,ld)); |
| 738 |
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.
|
31 | PetscCall(DSReallocateMat_Private(ds,DS_MAT_V,ld)); |
| 739 |
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.
|
31 | PetscCall(DSReallocateMat_Private(ds,DS_MAT_T,ld)); |
| 740 | |||
| 741 |
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.
|
31 | PetscCall(PetscMalloc1(ld,&ds->perm)); |
| 742 |
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.
|
31 | PetscCall(PetscArraycpy(ds->perm,perm,ds->ld)); |
| 743 |
6/8✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
31 | PetscCall(PetscFree(perm)); |
| 744 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
7 | PetscFunctionReturn(PETSC_SUCCESS); |
| 745 | } | ||
| 746 | |||
| 747 | /*MC | ||
| 748 | DSSVD - Dense Singular Value Decomposition. | ||
| 749 | |||
| 750 | Notes: | ||
| 751 | The problem is expressed as $A = U\Sigma V^*$, where $A$ is rectangular in | ||
| 752 | general, with $n$ rows and $m$ columns. $\Sigma$ is a diagonal matrix whose diagonal | ||
| 753 | elements are the arguments of `DSSolve()`. After solve, $A$ is overwritten | ||
| 754 | with $\Sigma$. | ||
| 755 | |||
| 756 | The orthogonal (or unitary) matrices of left and right singular vectors, $U$ | ||
| 757 | and $V$, have size $n$ and $m$, respectively. The number of columns $m$ must | ||
| 758 | be specified via `DSSVDSetDimensions()`. | ||
| 759 | |||
| 760 | If the `DS` object is in the intermediate state, $A$ is assumed to be in upper | ||
| 761 | bidiagonal form (possibly with an arrow) and is stored in compact format | ||
| 762 | on matrix $T$. Otherwise, no particular structure is assumed. The compact | ||
| 763 | storage is implemented for the square case only, $m=n$. The extra row should | ||
| 764 | be interpreted in this case as an extra column. | ||
| 765 | |||
| 766 | Used DS matrices: | ||
| 767 | + `DS_MAT_A` - problem matrix (used only if `compact=PETSC_FALSE`) | ||
| 768 | . `DS_MAT_T` - upper bidiagonal matrix | ||
| 769 | . `DS_MAT_U` - left singular vectors | ||
| 770 | - `DS_MAT_V` - right singular vectors | ||
| 771 | |||
| 772 | Implemented methods: | ||
| 773 | + 0 - Implicit zero-shift QR for bidiagonals (`_bdsqr`) | ||
| 774 | - 1 - Divide and Conquer (`_bdsdc` or `_gesdd`) | ||
| 775 | |||
| 776 | Level: beginner | ||
| 777 | |||
| 778 | .seealso: [](sec:ds), `DSCreate()`, `DSSetType()`, `DSType`, `DSSVDSetDimensions()`, `DSSetCompact()` | ||
| 779 | M*/ | ||
| 780 | 1001 | SLEPC_EXTERN PetscErrorCode DSCreate_SVD(DS ds) | |
| 781 | { | ||
| 782 | 1001 | DS_SVD *ctx; | |
| 783 | |||
| 784 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
1001 | PetscFunctionBegin; |
| 785 |
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.
|
1001 | PetscCall(PetscNew(&ctx)); |
| 786 | 1001 | ds->data = (void*)ctx; | |
| 787 | |||
| 788 | 1001 | ds->ops->allocate = DSAllocate_SVD; | |
| 789 | 1001 | ds->ops->view = DSView_SVD; | |
| 790 | 1001 | ds->ops->vectors = DSVectors_SVD; | |
| 791 | 1001 | ds->ops->solve[0] = DSSolve_SVD_QR; | |
| 792 | 1001 | ds->ops->solve[1] = DSSolve_SVD_DC; | |
| 793 | 1001 | ds->ops->sort = DSSort_SVD; | |
| 794 | 1001 | ds->ops->truncate = DSTruncate_SVD; | |
| 795 | 1001 | ds->ops->update = DSUpdateExtraRow_SVD; | |
| 796 | 1001 | ds->ops->destroy = DSDestroy_SVD; | |
| 797 | 1001 | ds->ops->matgetsize = DSMatGetSize_SVD; | |
| 798 | #if !defined(PETSC_HAVE_MPIUNI) | ||
| 799 | 1001 | ds->ops->synchronize = DSSynchronize_SVD; | |
| 800 | #endif | ||
| 801 | 1001 | ds->ops->setcompact = DSSetCompact_SVD; | |
| 802 | 1001 | ds->ops->reallocate = DSReallocate_SVD; | |
| 803 |
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.
|
1001 | PetscCall(PetscObjectComposeFunction((PetscObject)ds,"DSSVDSetDimensions_C",DSSVDSetDimensions_SVD)); |
| 804 |
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.
|
1001 | PetscCall(PetscObjectComposeFunction((PetscObject)ds,"DSSVDGetDimensions_C",DSSVDGetDimensions_SVD)); |
| 805 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
187 | PetscFunctionReturn(PETSC_SUCCESS); |
| 806 | } | ||
| 807 |