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> | ||
12 | #include <slepcblaslapack.h> | ||
13 | |||
14 | typedef struct { | ||
15 | PetscScalar *wr,*wi; /* eigenvalues of B */ | ||
16 | } DS_NHEPTS; | ||
17 | |||
18 | 159 | static PetscErrorCode DSAllocate_NHEPTS(DS ds,PetscInt ld) | |
19 | { | ||
20 | 159 | DS_NHEPTS *ctx = (DS_NHEPTS*)ds->data; | |
21 | |||
22 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
159 | PetscFunctionBegin; |
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.
|
159 | PetscCall(DSAllocateMat_Private(ds,DS_MAT_A)); |
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.
|
159 | PetscCall(DSAllocateMat_Private(ds,DS_MAT_B)); |
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.
|
159 | PetscCall(DSAllocateMat_Private(ds,DS_MAT_Q)); |
26 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
159 | PetscCall(DSAllocateMat_Private(ds,DS_MAT_Z)); |
27 |
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.
|
159 | PetscCall(PetscFree(ds->perm)); |
28 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
159 | PetscCall(PetscMalloc1(ld,&ds->perm)); |
29 |
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.
|
159 | PetscCall(PetscMalloc1(ld,&ctx->wr)); |
30 | #if !defined(PETSC_USE_COMPLEX) | ||
31 |
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.
|
89 | PetscCall(PetscMalloc1(ld,&ctx->wi)); |
32 | #endif | ||
33 |
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.
|
31 | PetscFunctionReturn(PETSC_SUCCESS); |
34 | } | ||
35 | |||
36 | ✗ | static PetscErrorCode DSView_NHEPTS(DS ds,PetscViewer viewer) | |
37 | { | ||
38 | ✗ | PetscViewerFormat format; | |
39 | |||
40 | ✗ | PetscFunctionBegin; | |
41 | ✗ | PetscCall(PetscViewerGetFormat(viewer,&format)); | |
42 | ✗ | if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) PetscFunctionReturn(PETSC_SUCCESS); | |
43 | ✗ | PetscCall(DSViewMat(ds,viewer,DS_MAT_A)); | |
44 | ✗ | PetscCall(DSViewMat(ds,viewer,DS_MAT_B)); | |
45 | ✗ | if (ds->state>DS_STATE_INTERMEDIATE) { | |
46 | ✗ | PetscCall(DSViewMat(ds,viewer,DS_MAT_Q)); | |
47 | ✗ | PetscCall(DSViewMat(ds,viewer,DS_MAT_Z)); | |
48 | } | ||
49 | ✗ | if (ds->omat[DS_MAT_X]) PetscCall(DSViewMat(ds,viewer,DS_MAT_X)); | |
50 | ✗ | if (ds->omat[DS_MAT_Y]) PetscCall(DSViewMat(ds,viewer,DS_MAT_Y)); | |
51 | ✗ | PetscFunctionReturn(PETSC_SUCCESS); | |
52 | } | ||
53 | |||
54 | 3884 | static PetscErrorCode DSVectors_NHEPTS_Eigen_Some(DS ds,PetscInt *k,PetscReal *rnorm,PetscBool left) | |
55 | { | ||
56 | 3884 | PetscInt i; | |
57 | 3884 | PetscBLASInt mm=1,mout,info,ld,n,*select,inc=1,cols=1,zero=0; | |
58 | 3884 | PetscScalar sone=1.0,szero=0.0; | |
59 | 3884 | PetscReal norm,done=1.0; | |
60 | 3884 | PetscBool iscomplex = PETSC_FALSE; | |
61 | 3884 | PetscScalar *X,*Y; | |
62 | 3884 | const PetscScalar *A,*Q; | |
63 | |||
64 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
3884 | PetscFunctionBegin; |
65 |
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.
|
3884 | PetscCall(PetscBLASIntCast(ds->n,&n)); |
66 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
3884 | PetscCall(PetscBLASIntCast(ds->ld,&ld)); |
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.
|
3884 | PetscCall(DSAllocateWork_Private(ds,0,0,ld)); |
68 | 3884 | select = ds->iwork; | |
69 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
81104 | for (i=0;i<n;i++) select[i] = (PetscBLASInt)PETSC_FALSE; |
70 | |||
71 | /* compute k-th eigenvector Y of A */ | ||
72 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
3884 | PetscCall(MatDenseGetArrayRead(ds->omat[left?DS_MAT_B:DS_MAT_A],&A)); |
73 |
7/8✓ 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 taken 2 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
5826 | PetscCall(MatDenseGetArray(ds->omat[left?DS_MAT_Y:DS_MAT_X],&X)); |
74 | 3884 | Y = X+(*k)*ld; | |
75 | 3884 | select[*k] = (PetscBLASInt)PETSC_TRUE; | |
76 | #if !defined(PETSC_USE_COMPLEX) | ||
77 |
3/4✓ Branch 0 taken 5 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 5 times.
✓ Branch 3 taken 5 times.
|
2088 | if ((*k)<n-1 && A[(*k)+1+(*k)*ld]!=0.0) iscomplex = PETSC_TRUE; |
78 | 2088 | mm = iscomplex? 2: 1; | |
79 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
2088 | if (iscomplex) select[(*k)+1] = (PetscBLASInt)PETSC_TRUE; |
80 |
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.
|
2088 | PetscCall(DSAllocateWork_Private(ds,3*ld,0,0)); |
81 |
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.
|
2088 | PetscCallBLAS("LAPACKtrevc",LAPACKtrevc_("R","S",select,&n,(PetscScalar*)A,&ld,Y,&ld,Y,&ld,&mm,&mout,ds->work,&info)); |
82 | #else | ||
83 |
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.
|
1796 | PetscCall(DSAllocateWork_Private(ds,2*ld,ld,0)); |
84 |
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.
|
1796 | PetscCallBLAS("LAPACKtrevc",LAPACKtrevc_("R","S",select,&n,(PetscScalar*)A,&ld,Y,&ld,Y,&ld,&mm,&mout,ds->work,ds->rwork,&info)); |
85 | #endif | ||
86 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
3884 | SlepcCheckLapackInfo("trevc",info); |
87 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
3884 | PetscCheck(mout==mm,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Inconsistent arguments"); |
88 |
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.
|
3884 | PetscCall(MatDenseRestoreArrayRead(ds->omat[left?DS_MAT_B:DS_MAT_A],&A)); |
89 | |||
90 | /* accumulate and normalize eigenvectors */ | ||
91 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
3884 | if (ds->state>=DS_STATE_CONDENSED) { |
92 |
7/8✓ 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 taken 2 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
5826 | PetscCall(MatDenseGetArrayRead(ds->omat[left?DS_MAT_Z:DS_MAT_Q],&Q)); |
93 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
3884 | PetscCall(PetscArraycpy(ds->work,Y,mout*ld)); |
94 |
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.
|
3884 | PetscCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&sone,Q,&ld,ds->work,&inc,&szero,Y,&inc)); |
95 | #if !defined(PETSC_USE_COMPLEX) | ||
96 |
12/22✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 1 times.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 1 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 1 times.
✓ Branch 10 taken 1 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 1 times.
✓ Branch 14 taken 1 times.
✗ Branch 15 not taken.
✓ Branch 16 taken 1 times.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
✓ Branch 19 taken 1 times.
✗ Branch 20 not taken.
✗ Branch 21 not taken.
|
2088 | if (iscomplex) PetscCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&sone,Q,&ld,ds->work+ld,&inc,&szero,Y+ld,&inc)); |
97 | #endif | ||
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.
|
3884 | PetscCall(MatDenseRestoreArrayRead(ds->omat[left?DS_MAT_Z:DS_MAT_Q],&Q)); |
99 | 3884 | cols = 1; | |
100 | 3884 | norm = BLASnrm2_(&n,Y,&inc); | |
101 | #if !defined(PETSC_USE_COMPLEX) | ||
102 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
2088 | if (iscomplex) { |
103 | 346 | norm = SlepcAbsEigenvalue(norm,BLASnrm2_(&n,Y+ld,&inc)); | |
104 | 346 | cols = 2; | |
105 | } | ||
106 | #endif | ||
107 |
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.
|
3884 | PetscCallBLAS("LAPACKlascl",LAPACKlascl_("G",&zero,&zero,&norm,&done,&n,&cols,Y,&ld,&info)); |
108 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
3884 | SlepcCheckLapackInfo("lascl",info); |
109 | } | ||
110 | |||
111 | /* set output arguments */ | ||
112 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
3884 | if (iscomplex) (*k)++; |
113 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
3884 | if (rnorm) { |
114 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
3884 | if (iscomplex) *rnorm = SlepcAbsEigenvalue(Y[n-1],Y[n-1+ld]); |
115 | 3538 | else *rnorm = PetscAbsScalar(Y[n-1]); | |
116 | } | ||
117 |
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.
|
3884 | PetscCall(MatDenseRestoreArray(ds->omat[left?DS_MAT_Y:DS_MAT_X],&X)); |
118 |
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.
|
764 | PetscFunctionReturn(PETSC_SUCCESS); |
119 | } | ||
120 | |||
121 | 318 | static PetscErrorCode DSVectors_NHEPTS_Eigen_All(DS ds,PetscBool left) | |
122 | { | ||
123 | 318 | PetscInt i; | |
124 | 318 | PetscBLASInt n,ld,mout,info,inc=1,cols,zero=0; | |
125 | 318 | PetscBool iscomplex; | |
126 | 318 | PetscScalar *X; | |
127 | 318 | const PetscScalar *A; | |
128 | 318 | PetscReal norm,done=1.0; | |
129 | 318 | const char *back; | |
130 | |||
131 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
318 | PetscFunctionBegin; |
132 |
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.
|
318 | PetscCall(PetscBLASIntCast(ds->n,&n)); |
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.
|
318 | PetscCall(PetscBLASIntCast(ds->ld,&ld)); |
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.
|
318 | PetscCall(MatDenseGetArrayRead(ds->omat[left?DS_MAT_B:DS_MAT_A],&A)); |
135 |
7/8✓ 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 taken 2 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
477 | PetscCall(MatDenseGetArrayWrite(ds->omat[left?DS_MAT_Y:DS_MAT_X],&X)); |
136 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
318 | if (ds->state>=DS_STATE_CONDENSED) { |
137 | /* DSSolve() has been called, backtransform with matrix Q */ | ||
138 | ✗ | back = "B"; | |
139 | ✗ | PetscCall(MatCopy(ds->omat[left?DS_MAT_Z:DS_MAT_Q],ds->omat[left?DS_MAT_Y:DS_MAT_X],SAME_NONZERO_PATTERN)); | |
140 | } else back = "A"; | ||
141 | #if !defined(PETSC_USE_COMPLEX) | ||
142 |
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.
|
178 | PetscCall(DSAllocateWork_Private(ds,3*ld,0,0)); |
143 |
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.
|
178 | PetscCallBLAS("LAPACKtrevc",LAPACKtrevc_("R",back,NULL,&n,(PetscScalar*)A,&ld,X,&ld,X,&ld,&n,&mout,ds->work,&info)); |
144 | #else | ||
145 |
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.
|
140 | PetscCall(DSAllocateWork_Private(ds,2*ld,ld,0)); |
146 |
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.
|
140 | PetscCallBLAS("LAPACKtrevc",LAPACKtrevc_("R",back,NULL,&n,(PetscScalar*)A,&ld,X,&ld,X,&ld,&n,&mout,ds->work,ds->rwork,&info)); |
147 | #endif | ||
148 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
318 | SlepcCheckLapackInfo("trevc",info); |
149 | |||
150 | /* normalize eigenvectors */ | ||
151 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1890 | for (i=0;i<n;i++) { |
152 |
4/4✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 5 times.
|
1572 | iscomplex = (i<n-1 && A[i+1+i*ld]!=0.0)? PETSC_TRUE: PETSC_FALSE; |
153 | 1572 | cols = 1; | |
154 | 1572 | norm = BLASnrm2_(&n,X+i*ld,&inc); | |
155 | #if !defined(PETSC_USE_COMPLEX) | ||
156 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
886 | if (iscomplex) { |
157 | 48 | norm = SlepcAbsEigenvalue(norm,BLASnrm2_(&n,X+(i+1)*ld,&inc)); | |
158 | 48 | cols = 2; | |
159 | } | ||
160 | #endif | ||
161 |
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.
|
1572 | PetscCallBLAS("LAPACKlascl",LAPACKlascl_("G",&zero,&zero,&norm,&done,&n,&cols,X+i*ld,&ld,&info)); |
162 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
1572 | SlepcCheckLapackInfo("lascl",info); |
163 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 10 times.
|
1572 | if (iscomplex) i++; |
164 | } | ||
165 |
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.
|
318 | PetscCall(MatDenseRestoreArrayRead(ds->omat[left?DS_MAT_B:DS_MAT_A],&A)); |
166 |
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.
|
318 | PetscCall(MatDenseRestoreArrayWrite(ds->omat[left?DS_MAT_Y:DS_MAT_X],&X)); |
167 |
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.
|
62 | PetscFunctionReturn(PETSC_SUCCESS); |
168 | } | ||
169 | |||
170 | 4202 | static PetscErrorCode DSVectors_NHEPTS(DS ds,DSMatType mat,PetscInt *j,PetscReal *rnorm) | |
171 | { | ||
172 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
4202 | PetscFunctionBegin; |
173 |
2/4✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
4202 | switch (mat) { |
174 | 2101 | case DS_MAT_X: | |
175 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
2101 | PetscCheck(!ds->refined,PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented yet"); |
176 |
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.
|
2101 | if (j) PetscCall(DSVectors_NHEPTS_Eigen_Some(ds,j,rnorm,PETSC_FALSE)); |
177 |
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.
|
159 | else PetscCall(DSVectors_NHEPTS_Eigen_All(ds,PETSC_FALSE)); |
178 | break; | ||
179 | 2101 | case DS_MAT_Y: | |
180 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
2101 | PetscCheck(!ds->refined,PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented yet"); |
181 |
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.
|
2101 | if (j) PetscCall(DSVectors_NHEPTS_Eigen_Some(ds,j,rnorm,PETSC_TRUE)); |
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.
|
159 | else PetscCall(DSVectors_NHEPTS_Eigen_All(ds,PETSC_TRUE)); |
183 | break; | ||
184 | ✗ | case DS_MAT_U: | |
185 | case DS_MAT_V: | ||
186 | ✗ | SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented yet"); | |
187 | ✗ | default: | |
188 | ✗ | SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter"); | |
189 | } | ||
190 |
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.
|
826 | PetscFunctionReturn(PETSC_SUCCESS); |
191 | } | ||
192 | |||
193 | 1006 | static PetscErrorCode DSSort_NHEPTS(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k) | |
194 | { | ||
195 | 1006 | DS_NHEPTS *ctx = (DS_NHEPTS*)ds->data; | |
196 | 1006 | PetscInt i,j,cont,id=0,*p,*idx,*idx2; | |
197 | 1006 | PetscReal s,t; | |
198 | #if defined(PETSC_USE_COMPLEX) | ||
199 | 475 | Mat A,U; | |
200 | #endif | ||
201 | |||
202 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
1006 | PetscFunctionBegin; |
203 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
1006 | PetscCheck(!rr || wr==rr,PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented yet"); |
204 |
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.
|
1006 | PetscCall(PetscMalloc3(ds->ld,&idx,ds->ld,&idx2,ds->ld,&p)); |
205 |
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.
|
1006 | PetscCall(DSSort_NHEP_Total(ds,DS_MAT_A,DS_MAT_Q,wr,wi)); |
206 | #if defined(PETSC_USE_COMPLEX) | ||
207 |
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.
|
475 | PetscCall(DSGetMat(ds,DS_MAT_B,&A)); |
208 |
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.
|
475 | PetscCall(MatConjugate(A)); |
209 |
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.
|
475 | PetscCall(DSRestoreMat(ds,DS_MAT_B,&A)); |
210 |
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.
|
475 | PetscCall(DSGetMat(ds,DS_MAT_Z,&U)); |
211 |
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.
|
475 | PetscCall(MatConjugate(U)); |
212 |
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.
|
475 | PetscCall(DSRestoreMat(ds,DS_MAT_Z,&U)); |
213 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
9538 | for (i=0;i<ds->n;i++) ctx->wr[i] = PetscConj(ctx->wr[i]); |
214 | #endif | ||
215 |
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.
|
1006 | PetscCall(DSSort_NHEP_Total(ds,DS_MAT_B,DS_MAT_Z,ctx->wr,ctx->wi)); |
216 | /* check correct eigenvalue correspondence */ | ||
217 | cont = 0; | ||
218 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
20702 | for (i=0;i<ds->n;i++) { |
219 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
19696 | if (SlepcAbsEigenvalue(ctx->wr[i]-wr[i],ctx->wi[i]-wi[i])>PETSC_SQRT_MACHINE_EPSILON) {idx2[cont] = i; idx[cont++] = i;} |
220 | 19696 | p[i] = -1; | |
221 | } | ||
222 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1006 | if (cont) { |
223 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
3502 | for (i=0;i<cont;i++) { |
224 | t = PETSC_MAX_REAL; | ||
225 |
6/6✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
✓ Branch 4 taken 10 times.
✓ Branch 5 taken 10 times.
|
85314 | for (j=0;j<cont;j++) if (idx2[j]!=-1 && (s=SlepcAbsEigenvalue(ctx->wr[idx[j]]-wr[idx[i]],ctx->wi[idx[j]]-wi[idx[i]]))<t) { id = j; t = s; } |
226 | 3325 | p[idx[i]] = idx[id]; | |
227 | 3325 | idx2[id] = -1; | |
228 | } | ||
229 |
4/4✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
|
4472 | for (i=0;i<ds->n;i++) if (p[i]==-1) p[i] = i; |
230 |
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.
|
177 | PetscCall(DSSortWithPermutation_NHEP_Private(ds,p,DS_MAT_B,DS_MAT_Z,ctx->wr,ctx->wi)); |
231 | } | ||
232 | #if defined(PETSC_USE_COMPLEX) | ||
233 |
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.
|
475 | PetscCall(DSGetMat(ds,DS_MAT_B,&A)); |
234 |
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.
|
475 | PetscCall(MatConjugate(A)); |
235 |
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.
|
475 | PetscCall(DSRestoreMat(ds,DS_MAT_B,&A)); |
236 |
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.
|
475 | PetscCall(DSGetMat(ds,DS_MAT_Z,&U)); |
237 |
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.
|
475 | PetscCall(MatConjugate(U)); |
238 |
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.
|
475 | PetscCall(DSRestoreMat(ds,DS_MAT_Z,&U)); |
239 | #endif | ||
240 |
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.
|
1006 | PetscCall(PetscFree3(idx,idx2,p)); |
241 |
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.
|
200 | PetscFunctionReturn(PETSC_SUCCESS); |
242 | } | ||
243 | |||
244 | 1006 | static PetscErrorCode DSUpdateExtraRow_NHEPTS(DS ds) | |
245 | { | ||
246 | 1006 | PetscInt i; | |
247 | 1006 | PetscBLASInt n,ld,incx=1; | |
248 | 1006 | PetscScalar *A,*x,*y,one=1.0,zero=0.0; | |
249 | 1006 | const PetscScalar *Q; | |
250 | |||
251 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
1006 | PetscFunctionBegin; |
252 |
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.
|
1006 | PetscCall(PetscBLASIntCast(ds->n,&n)); |
253 |
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.
|
1006 | PetscCall(PetscBLASIntCast(ds->ld,&ld)); |
254 |
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.
|
1006 | PetscCall(DSAllocateWork_Private(ds,2*ld,0,0)); |
255 | 1006 | x = ds->work; | |
256 | 1006 | y = ds->work+ld; | |
257 |
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.
|
1006 | PetscCall(MatDenseGetArray(ds->omat[DS_MAT_A],&A)); |
258 |
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.
|
1006 | PetscCall(MatDenseGetArrayRead(ds->omat[DS_MAT_Q],&Q)); |
259 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
20702 | for (i=0;i<n;i++) x[i] = PetscConj(A[n+i*ld]); |
260 |
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.
|
1006 | PetscCallBLAS("BLASgemv",BLASgemv_("C",&n,&n,&one,Q,&ld,x,&incx,&zero,y,&incx)); |
261 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
20702 | for (i=0;i<n;i++) A[n+i*ld] = PetscConj(y[i]); |
262 |
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.
|
1006 | PetscCall(MatDenseRestoreArray(ds->omat[DS_MAT_A],&A)); |
263 |
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.
|
1006 | PetscCall(MatDenseRestoreArrayRead(ds->omat[DS_MAT_Q],&Q)); |
264 |
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.
|
1006 | PetscCall(MatDenseGetArray(ds->omat[DS_MAT_B],&A)); |
265 |
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.
|
1006 | PetscCall(MatDenseGetArrayRead(ds->omat[DS_MAT_Z],&Q)); |
266 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
20702 | for (i=0;i<n;i++) x[i] = PetscConj(A[n+i*ld]); |
267 |
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.
|
1006 | PetscCallBLAS("BLASgemv",BLASgemv_("C",&n,&n,&one,Q,&ld,x,&incx,&zero,y,&incx)); |
268 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
20702 | for (i=0;i<n;i++) A[n+i*ld] = PetscConj(y[i]); |
269 |
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.
|
1006 | PetscCall(MatDenseRestoreArray(ds->omat[DS_MAT_B],&A)); |
270 |
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.
|
1006 | PetscCall(MatDenseRestoreArrayRead(ds->omat[DS_MAT_Z],&Q)); |
271 | 1006 | ds->k = n; | |
272 |
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.
|
1006 | PetscFunctionReturn(PETSC_SUCCESS); |
273 | } | ||
274 | |||
275 | 1006 | static PetscErrorCode DSSolve_NHEPTS(DS ds,PetscScalar *wr,PetscScalar *wi) | |
276 | { | ||
277 | 1006 | DS_NHEPTS *ctx = (DS_NHEPTS*)ds->data; | |
278 | |||
279 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
1006 | PetscFunctionBegin; |
280 | #if !defined(PETSC_USE_COMPLEX) | ||
281 |
2/8✗ Branch 0 not taken.
✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
531 | PetscAssertPointer(wi,3); |
282 | #endif | ||
283 |
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.
|
1006 | PetscCall(DSSolve_NHEP_Private(ds,DS_MAT_A,DS_MAT_Q,wr,wi)); |
284 |
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.
|
1006 | PetscCall(DSSolve_NHEP_Private(ds,DS_MAT_B,DS_MAT_Z,ctx->wr,ctx->wi)); |
285 |
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.
|
200 | PetscFunctionReturn(PETSC_SUCCESS); |
286 | } | ||
287 | |||
288 | #if !defined(PETSC_HAVE_MPIUNI) | ||
289 | 280 | static PetscErrorCode DSSynchronize_NHEPTS(DS ds,PetscScalar eigr[],PetscScalar eigi[]) | |
290 | { | ||
291 | 280 | PetscInt ld=ds->ld,l=ds->l,k; | |
292 | 280 | PetscMPIInt n,rank,off=0,size,ldn; | |
293 | 280 | DS_NHEPTS *ctx = (DS_NHEPTS*)ds->data; | |
294 | 280 | PetscScalar *A,*B,*Q,*Z; | |
295 | |||
296 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
280 | PetscFunctionBegin; |
297 | 280 | k = 2*(ds->n-l)*ld; | |
298 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
280 | if (ds->state>DS_STATE_RAW) k += 2*(ds->n-l)*ld; |
299 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
280 | if (eigr) k += ds->n-l; |
300 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
280 | if (eigi) k += ds->n-l; |
301 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
280 | if (ctx->wr) k += ds->n-l; |
302 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
280 | if (ctx->wi) k += ds->n-l; |
303 |
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.
|
280 | PetscCall(DSAllocateWork_Private(ds,k,0,0)); |
304 |
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.
|
280 | PetscCall(PetscMPIIntCast(k*sizeof(PetscScalar),&size)); |
305 |
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.
|
280 | PetscCall(PetscMPIIntCast(ds->n-l,&n)); |
306 |
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.
|
280 | PetscCall(PetscMPIIntCast(ld*(ds->n-l),&ldn)); |
307 |
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.
|
280 | PetscCall(MatDenseGetArray(ds->omat[DS_MAT_A],&A)); |
308 |
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.
|
280 | PetscCall(MatDenseGetArray(ds->omat[DS_MAT_B],&B)); |
309 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
280 | if (ds->state>DS_STATE_RAW) { |
310 |
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.
|
280 | PetscCall(MatDenseGetArray(ds->omat[DS_MAT_Q],&Q)); |
311 |
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.
|
280 | PetscCall(MatDenseGetArray(ds->omat[DS_MAT_Z],&Z)); |
312 | } | ||
313 |
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.
|
280 | PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)ds),&rank)); |
314 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
280 | if (!rank) { |
315 |
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.
|
140 | PetscCallMPI(MPI_Pack(A+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds))); |
316 |
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.
|
140 | PetscCallMPI(MPI_Pack(B+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds))); |
317 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
140 | if (ds->state>DS_STATE_RAW) { |
318 |
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.
|
140 | PetscCallMPI(MPI_Pack(Q+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds))); |
319 |
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.
|
140 | PetscCallMPI(MPI_Pack(Z+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds))); |
320 | } | ||
321 |
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.
|
140 | if (eigr) PetscCallMPI(MPI_Pack(eigr+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds))); |
322 | #if !defined(PETSC_USE_COMPLEX) | ||
323 |
15/30✓ Branch 0 taken 5 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 1 times.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 1 times.
✓ Branch 5 taken 4 times.
✓ Branch 6 taken 1 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 1 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 1 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 1 times.
✓ Branch 14 taken 1 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 1 times.
✓ Branch 18 taken 1 times.
✗ Branch 19 not taken.
✓ Branch 20 taken 1 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 1 times.
✗ Branch 24 not taken.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 1 times.
✗ Branch 28 not taken.
✗ Branch 29 not taken.
|
70 | if (eigi) PetscCallMPI(MPI_Pack(eigi+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds))); |
324 | #endif | ||
325 |
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.
|
140 | if (ctx->wr) PetscCallMPI(MPI_Pack(ctx->wr+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds))); |
326 |
16/30✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 1 times.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 1 times.
✓ Branch 5 taken 4 times.
✓ Branch 6 taken 1 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 1 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 1 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 1 times.
✓ Branch 14 taken 1 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 1 times.
✓ Branch 18 taken 1 times.
✗ Branch 19 not taken.
✓ Branch 20 taken 1 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 1 times.
✗ Branch 24 not taken.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 1 times.
✗ Branch 28 not taken.
✗ Branch 29 not taken.
|
140 | if (ctx->wi) PetscCallMPI(MPI_Pack(ctx->wi+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds))); |
327 | } | ||
328 |
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.
|
560 | PetscCallMPI(MPI_Bcast(ds->work,size,MPI_BYTE,0,PetscObjectComm((PetscObject)ds))); |
329 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
280 | if (rank) { |
330 |
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.
|
140 | PetscCallMPI(MPI_Unpack(ds->work,size,&off,A+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds))); |
331 |
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.
|
140 | PetscCallMPI(MPI_Unpack(ds->work,size,&off,B+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds))); |
332 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
140 | if (ds->state>DS_STATE_RAW) { |
333 |
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.
|
140 | PetscCallMPI(MPI_Unpack(ds->work,size,&off,Q+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds))); |
334 |
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.
|
140 | PetscCallMPI(MPI_Unpack(ds->work,size,&off,Z+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds))); |
335 | } | ||
336 |
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.
|
140 | if (eigr) PetscCallMPI(MPI_Unpack(ds->work,size,&off,eigr+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds))); |
337 | #if !defined(PETSC_USE_COMPLEX) | ||
338 |
15/30✓ Branch 0 taken 5 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 1 times.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 1 times.
✓ Branch 5 taken 4 times.
✓ Branch 6 taken 1 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 1 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 1 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 1 times.
✓ Branch 14 taken 1 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 1 times.
✓ Branch 18 taken 1 times.
✗ Branch 19 not taken.
✓ Branch 20 taken 1 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 1 times.
✗ Branch 24 not taken.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 1 times.
✗ Branch 28 not taken.
✗ Branch 29 not taken.
|
70 | if (eigi) PetscCallMPI(MPI_Unpack(ds->work,size,&off,eigi+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds))); |
339 | #endif | ||
340 |
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.
|
140 | if (ctx->wr) PetscCallMPI(MPI_Unpack(ds->work,size,&off,ctx->wr+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds))); |
341 |
16/30✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 1 times.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 1 times.
✓ Branch 5 taken 4 times.
✓ Branch 6 taken 1 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 1 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 1 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 1 times.
✓ Branch 14 taken 1 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 1 times.
✓ Branch 18 taken 1 times.
✗ Branch 19 not taken.
✓ Branch 20 taken 1 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 1 times.
✗ Branch 24 not taken.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 1 times.
✗ Branch 28 not taken.
✗ Branch 29 not taken.
|
140 | if (ctx->wi) PetscCallMPI(MPI_Unpack(ds->work,size,&off,ctx->wi+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds))); |
342 | } | ||
343 |
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.
|
280 | PetscCall(MatDenseRestoreArray(ds->omat[DS_MAT_A],&A)); |
344 |
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.
|
280 | PetscCall(MatDenseRestoreArray(ds->omat[DS_MAT_B],&B)); |
345 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
280 | if (ds->state>DS_STATE_RAW) { |
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.
|
280 | PetscCall(MatDenseRestoreArray(ds->omat[DS_MAT_Q],&Q)); |
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.
|
280 | PetscCall(MatDenseRestoreArray(ds->omat[DS_MAT_Z],&Z)); |
348 | } | ||
349 |
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.
|
56 | PetscFunctionReturn(PETSC_SUCCESS); |
350 | } | ||
351 | #endif | ||
352 | |||
353 | 847 | static PetscErrorCode DSGetTruncateSize_NHEPTS(DS ds,PetscInt l,PetscInt n,PetscInt *k) | |
354 | { | ||
355 | #if !defined(PETSC_USE_COMPLEX) | ||
356 | 442 | const PetscScalar *A,*B; | |
357 | #endif | ||
358 | |||
359 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
847 | PetscFunctionBegin; |
360 | #if !defined(PETSC_USE_COMPLEX) | ||
361 |
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.
|
442 | PetscCall(MatDenseGetArrayRead(ds->omat[DS_MAT_A],&A)); |
362 |
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.
|
442 | PetscCall(MatDenseGetArrayRead(ds->omat[DS_MAT_B],&B)); |
363 |
3/4✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 5 times.
|
442 | if (A[l+(*k)+(l+(*k)-1)*ds->ld] != 0.0 || B[l+(*k)+(l+(*k)-1)*ds->ld] != 0.0) { |
364 |
1/2✓ Branch 0 taken 5 times.
✗ Branch 1 not taken.
|
104 | if (l+(*k)<n-1) (*k)++; |
365 | ✗ | else (*k)--; | |
366 | } | ||
367 |
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.
|
442 | PetscCall(MatDenseRestoreArrayRead(ds->omat[DS_MAT_A],&A)); |
368 |
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.
|
442 | PetscCall(MatDenseRestoreArrayRead(ds->omat[DS_MAT_B],&B)); |
369 | #endif | ||
370 |
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.
|
483 | PetscFunctionReturn(PETSC_SUCCESS); |
371 | } | ||
372 | |||
373 | 1006 | static PetscErrorCode DSTruncate_NHEPTS(DS ds,PetscInt n,PetscBool trim) | |
374 | { | ||
375 | 1006 | PetscInt i,ld=ds->ld,l=ds->l; | |
376 | 1006 | PetscScalar *A,*B; | |
377 | |||
378 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
1006 | PetscFunctionBegin; |
379 |
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.
|
1006 | PetscCall(MatDenseGetArray(ds->omat[DS_MAT_A],&A)); |
380 |
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.
|
1006 | PetscCall(MatDenseGetArray(ds->omat[DS_MAT_B],&B)); |
381 | #if defined(PETSC_USE_DEBUG) | ||
382 | /* make sure diagonal 2x2 block is not broken */ | ||
383 |
4/12✓ 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 not taken.
✓ Branch 7 taken 2 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
|
200 | PetscCheck(ds->state<DS_STATE_CONDENSED || n==0 || n==ds->n || A[n+(n-1)*ld]==0.0 || B[n+(n-1)*ld]==0.0,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"The given size would break a 2x2 block, call DSGetTruncateSize() first"); |
384 | #endif | ||
385 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1006 | if (trim) { |
386 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
159 | if (ds->extrarow) { /* clean extra row */ |
387 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
3004 | for (i=l;i<ds->n;i++) { A[ds->n+i*ld] = 0.0; B[ds->n+i*ld] = 0.0; } |
388 | } | ||
389 | 159 | ds->l = 0; | |
390 | 159 | ds->k = 0; | |
391 | 159 | ds->n = n; | |
392 | 159 | ds->t = ds->n; /* truncated length equal to the new dimension */ | |
393 | } else { | ||
394 |
2/4✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
|
847 | if (ds->extrarow && ds->k==ds->n) { |
395 | /* copy entries of extra row to the new position, then clean last row */ | ||
396 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
8754 | for (i=l;i<n;i++) { A[n+i*ld] = A[ds->n+i*ld]; B[n+i*ld] = B[ds->n+i*ld]; } |
397 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
16900 | for (i=l;i<ds->n;i++) { A[ds->n+i*ld] = 0.0; B[ds->n+i*ld] = 0.0; } |
398 | } | ||
399 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
847 | ds->k = (ds->extrarow)? n: 0; |
400 | 847 | ds->t = ds->n; /* truncated length equal to previous dimension */ | |
401 | 847 | ds->n = n; | |
402 | } | ||
403 |
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.
|
1006 | PetscCall(MatDenseRestoreArray(ds->omat[DS_MAT_A],&A)); |
404 |
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.
|
1006 | PetscCall(MatDenseRestoreArray(ds->omat[DS_MAT_B],&B)); |
405 |
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.
|
200 | PetscFunctionReturn(PETSC_SUCCESS); |
406 | } | ||
407 | |||
408 | 159 | static PetscErrorCode DSDestroy_NHEPTS(DS ds) | |
409 | { | ||
410 | 159 | DS_NHEPTS *ctx = (DS_NHEPTS*)ds->data; | |
411 | |||
412 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
159 | PetscFunctionBegin; |
413 |
6/10✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✓ Branch 5 taken 8 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
|
159 | if (ctx->wr) PetscCall(PetscFree(ctx->wr)); |
414 |
7/10✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 5 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 1 times.
✓ Branch 5 taken 4 times.
✓ Branch 6 taken 1 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 1 times.
|
159 | if (ctx->wi) PetscCall(PetscFree(ctx->wi)); |
415 |
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.
|
159 | PetscCall(PetscFree(ds->data)); |
416 |
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.
|
31 | PetscFunctionReturn(PETSC_SUCCESS); |
417 | } | ||
418 | |||
419 | 6242 | static PetscErrorCode DSMatGetSize_NHEPTS(DS ds,DSMatType t,PetscInt *rows,PetscInt *cols) | |
420 | { | ||
421 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
6242 | PetscFunctionBegin; |
422 |
3/4✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
|
6242 | *rows = ((t==DS_MAT_A || t==DS_MAT_B) && ds->extrarow)? ds->n+1: ds->n; |
423 | 6242 | *cols = ds->n; | |
424 |
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.
|
6242 | PetscFunctionReturn(PETSC_SUCCESS); |
425 | } | ||
426 | |||
427 | /*MC | ||
428 | DSNHEPTS - Dense Non-Hermitian Eigenvalue Problem (special variant intended | ||
429 | for two-sided Krylov solvers). | ||
430 | |||
431 | Level: beginner | ||
432 | |||
433 | Notes: | ||
434 | Two related problems are solved, A*X = X*Lambda and B*Y = Y*Lambda', where A and | ||
435 | B are supposed to come from the Arnoldi factorizations of a certain matrix and its | ||
436 | (conjugate) transpose, respectively. Hence, in exact arithmetic the columns of Y | ||
437 | are equal to the left eigenvectors of A. Lambda is a diagonal matrix whose diagonal | ||
438 | elements are the arguments of DSSolve(). After solve, A is overwritten with the | ||
439 | upper quasi-triangular matrix T of the (real) Schur form, A*Q = Q*T, and similarly | ||
440 | another (real) Schur relation is computed, B*Z = Z*S, overwriting B. | ||
441 | |||
442 | In the intermediate state A and B are reduced to upper Hessenberg form. | ||
443 | |||
444 | When left eigenvectors DS_MAT_Y are requested, right eigenvectors of B are returned, | ||
445 | while DS_MAT_X contains right eigenvectors of A. | ||
446 | |||
447 | Used DS matrices: | ||
448 | + DS_MAT_A - first problem matrix obtained from Arnoldi | ||
449 | . DS_MAT_B - second problem matrix obtained from Arnoldi on the transpose | ||
450 | . DS_MAT_Q - orthogonal/unitary transformation that reduces A to Hessenberg form | ||
451 | (intermediate step) or matrix of orthogonal Schur vectors of A | ||
452 | - DS_MAT_Z - orthogonal/unitary transformation that reduces B to Hessenberg form | ||
453 | (intermediate step) or matrix of orthogonal Schur vectors of B | ||
454 | |||
455 | Implemented methods: | ||
456 | . 0 - Implicit QR (_hseqr) | ||
457 | |||
458 | .seealso: DSCreate(), DSSetType(), DSType | ||
459 | M*/ | ||
460 | 159 | SLEPC_EXTERN PetscErrorCode DSCreate_NHEPTS(DS ds) | |
461 | { | ||
462 | 159 | DS_NHEPTS *ctx; | |
463 | |||
464 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
159 | PetscFunctionBegin; |
465 |
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.
|
159 | PetscCall(PetscNew(&ctx)); |
466 | 159 | ds->data = (void*)ctx; | |
467 | |||
468 | 159 | ds->ops->allocate = DSAllocate_NHEPTS; | |
469 | 159 | ds->ops->view = DSView_NHEPTS; | |
470 | 159 | ds->ops->vectors = DSVectors_NHEPTS; | |
471 | 159 | ds->ops->solve[0] = DSSolve_NHEPTS; | |
472 | 159 | ds->ops->sort = DSSort_NHEPTS; | |
473 | #if !defined(PETSC_HAVE_MPIUNI) | ||
474 | 159 | ds->ops->synchronize = DSSynchronize_NHEPTS; | |
475 | #endif | ||
476 | 159 | ds->ops->gettruncatesize = DSGetTruncateSize_NHEPTS; | |
477 | 159 | ds->ops->truncate = DSTruncate_NHEPTS; | |
478 | 159 | ds->ops->update = DSUpdateExtraRow_NHEPTS; | |
479 | 159 | ds->ops->destroy = DSDestroy_NHEPTS; | |
480 | 159 | ds->ops->matgetsize = DSMatGetSize_NHEPTS; | |
481 |
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.
|
159 | PetscFunctionReturn(PETSC_SUCCESS); |
482 | } | ||
483 |