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 | This file implements a wrapper to eigensolvers in Elemental. | ||
12 | */ | ||
13 | |||
14 | #include <slepc/private/epsimpl.h> /*I "slepceps.h" I*/ | ||
15 | #include <petsc/private/petscelemental.h> | ||
16 | |||
17 | typedef struct { | ||
18 | Mat Ae,Be; /* converted matrices */ | ||
19 | } EPS_Elemental; | ||
20 | |||
21 | 19 | static PetscErrorCode EPSSetUp_Elemental(EPS eps) | |
22 | { | ||
23 | 19 | EPS_Elemental *ctx = (EPS_Elemental*)eps->data; | |
24 | 19 | Mat A,B; | |
25 | 19 | PetscInt nmat; | |
26 | 19 | PetscBool isshift; | |
27 | 19 | PetscScalar shift; | |
28 | |||
29 | 19 | PetscFunctionBegin; | |
30 |
4/10✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
|
19 | EPSCheckHermitianDefinite(eps); |
31 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
19 | EPSCheckNotStructured(eps); |
32 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
19 | PetscCall(PetscObjectTypeCompare((PetscObject)eps->st,STSHIFT,&isshift)); |
33 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
19 | PetscCheck(isshift,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support spectral transformations"); |
34 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
19 | if (eps->nev==0) eps->nev = 1; |
35 | 19 | eps->ncv = eps->n; | |
36 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
19 | if (eps->mpd!=PETSC_DETERMINE) PetscCall(PetscInfo(eps,"Warning: parameter mpd ignored\n")); |
37 |
2/2✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
|
19 | if (eps->max_it==PETSC_DETERMINE) eps->max_it = 1; |
38 |
3/4✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
|
19 | if (!eps->which) PetscCall(EPSSetWhichEigenpairs_Default(eps)); |
39 |
1/6✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
19 | PetscCheck(eps->which!=EPS_ALL || eps->inta==eps->intb,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support interval computation"); |
40 |
4/14✗ 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 taken 2 times.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
|
19 | EPSCheckUnsupported(eps,EPS_FEATURE_BALANCE | EPS_FEATURE_ARBITRARY | EPS_FEATURE_REGION); |
41 |
5/12✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
|
19 | EPSCheckIgnored(eps,EPS_FEATURE_EXTRACTION | EPS_FEATURE_CONVERGENCE | EPS_FEATURE_STOPPING); |
42 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
19 | PetscCall(EPSAllocateSolution(eps,0)); |
43 | |||
44 | /* convert matrices */ | ||
45 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
19 | PetscCall(MatDestroy(&ctx->Ae)); |
46 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
19 | PetscCall(MatDestroy(&ctx->Be)); |
47 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
19 | PetscCall(STGetNumMatrices(eps->st,&nmat)); |
48 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
19 | PetscCall(STGetMatrix(eps->st,0,&A)); |
49 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
19 | PetscCall(MatConvert(A,MATELEMENTAL,MAT_INITIAL_MATRIX,&ctx->Ae)); |
50 |
2/2✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
|
19 | if (nmat>1) { |
51 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
6 | PetscCall(STGetMatrix(eps->st,1,&B)); |
52 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
6 | PetscCall(MatConvert(B,MATELEMENTAL,MAT_INITIAL_MATRIX,&ctx->Be)); |
53 | } | ||
54 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
19 | PetscCall(STGetShift(eps->st,&shift)); |
55 |
2/2✓ Branch 0 taken 1 times.
✓ Branch 1 taken 1 times.
|
19 | if (shift != 0.0) { |
56 | ✗ | if (nmat>1) PetscCall(MatAXPY(ctx->Ae,-shift,ctx->Be,SAME_NONZERO_PATTERN)); | |
57 | ✗ | else PetscCall(MatShift(ctx->Ae,-shift)); | |
58 | } | ||
59 | PetscFunctionReturn(PETSC_SUCCESS); | ||
60 | } | ||
61 | |||
62 | 19 | static PetscErrorCode EPSSolve_Elemental(EPS eps) | |
63 | { | ||
64 | 19 | EPS_Elemental *ctx = (EPS_Elemental*)eps->data; | |
65 | 19 | Mat A = ctx->Ae,B = ctx->Be,Q,V; | |
66 | 19 | Mat_Elemental *a = (Mat_Elemental*)A->data,*b,*q; | |
67 | 19 | PetscInt i,rrank,ridx,erow; | |
68 | |||
69 | 19 | PetscFunctionBegin; | |
70 | 19 | El::DistMatrix<PetscReal,El::VR,El::STAR> w(*a->grid); | |
71 |
2/6✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
19 | PetscCall(MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&Q)); |
72 | 19 | q = (Mat_Elemental*)Q->data; | |
73 | |||
74 |
2/2✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
|
19 | if (B) { |
75 | 6 | b = (Mat_Elemental*)B->data; | |
76 |
2/4✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
|
6 | El::HermitianGenDefEig(El::AXBX,El::LOWER,*a->emat,*b->emat,w,*q->emat); |
77 |
2/4✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
|
13 | } else El::HermitianEig(El::LOWER,*a->emat,w,*q->emat); |
78 | |||
79 |
2/2✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
|
3405 | for (i=0;i<eps->ncv;i++) { |
80 | 3386 | P2RO(A,0,i,&rrank,&ridx); | |
81 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
3386 | RO2E(A,0,rrank,ridx,&erow); |
82 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
3386 | eps->eigr[i] = w.Get(erow,0); |
83 | } | ||
84 |
2/6✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
19 | PetscCall(BVGetMat(eps->V,&V)); |
85 |
2/6✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
19 | PetscCall(MatConvert(Q,MATDENSE,MAT_REUSE_MATRIX,&V)); |
86 |
2/6✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
19 | PetscCall(BVRestoreMat(eps->V,&V)); |
87 |
2/6✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
19 | PetscCall(MatDestroy(&Q)); |
88 | |||
89 | 19 | eps->nconv = eps->ncv; | |
90 | 19 | eps->its = 1; | |
91 | 19 | eps->reason = EPS_CONVERGED_TOL; | |
92 | 19 | PetscFunctionReturn(PETSC_SUCCESS); | |
93 | 19 | } | |
94 | |||
95 | 13 | static PetscErrorCode EPSDestroy_Elemental(EPS eps) | |
96 | { | ||
97 | 13 | PetscFunctionBegin; | |
98 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
13 | PetscCall(PetscFree(eps->data)); |
99 | PetscFunctionReturn(PETSC_SUCCESS); | ||
100 | } | ||
101 | |||
102 | 13 | static PetscErrorCode EPSReset_Elemental(EPS eps) | |
103 | { | ||
104 | 13 | EPS_Elemental *ctx = (EPS_Elemental*)eps->data; | |
105 | |||
106 | 13 | PetscFunctionBegin; | |
107 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
13 | PetscCall(MatDestroy(&ctx->Ae)); |
108 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
13 | PetscCall(MatDestroy(&ctx->Be)); |
109 | PetscFunctionReturn(PETSC_SUCCESS); | ||
110 | } | ||
111 | |||
112 | 13 | SLEPC_EXTERN PetscErrorCode EPSCreate_Elemental(EPS eps) | |
113 | { | ||
114 | 13 | EPS_Elemental *ctx; | |
115 | |||
116 | 13 | PetscFunctionBegin; | |
117 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
13 | PetscCall(PetscNew(&ctx)); |
118 | 13 | eps->data = (void*)ctx; | |
119 | |||
120 | 13 | eps->categ = EPS_CATEGORY_OTHER; | |
121 | |||
122 | 13 | eps->ops->solve = EPSSolve_Elemental; | |
123 | 13 | eps->ops->setup = EPSSetUp_Elemental; | |
124 | 13 | eps->ops->setupsort = EPSSetUpSort_Basic; | |
125 | 13 | eps->ops->destroy = EPSDestroy_Elemental; | |
126 | 13 | eps->ops->reset = EPSReset_Elemental; | |
127 | 13 | eps->ops->backtransform = EPSBackTransform_Default; | |
128 | 13 | eps->ops->setdefaultst = EPSSetDefaultST_NoFactor; | |
129 | 13 | PetscFunctionReturn(PETSC_SUCCESS); | |
130 | } | ||
131 |