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 | #pragma once | ||
12 | |||
13 | #include <slepcmfn.h> | ||
14 | #include <slepc/private/slepcimpl.h> | ||
15 | |||
16 | /* SUBMANSEC = MFN */ | ||
17 | |||
18 | SLEPC_EXTERN PetscBool MFNRegisterAllCalled; | ||
19 | SLEPC_EXTERN PetscBool MFNMonitorRegisterAllCalled; | ||
20 | SLEPC_EXTERN PetscErrorCode MFNRegisterAll(void); | ||
21 | SLEPC_EXTERN PetscErrorCode MFNMonitorRegisterAll(void); | ||
22 | SLEPC_EXTERN PetscLogEvent MFN_SetUp,MFN_Solve; | ||
23 | |||
24 | typedef struct _MFNOps *MFNOps; | ||
25 | |||
26 | struct _MFNOps { | ||
27 | PetscErrorCode (*solve)(MFN,Vec,Vec); | ||
28 | PetscErrorCode (*setup)(MFN); | ||
29 | PetscErrorCode (*setfromoptions)(MFN,PetscOptionItems); | ||
30 | PetscErrorCode (*publishoptions)(MFN); | ||
31 | PetscErrorCode (*destroy)(MFN); | ||
32 | PetscErrorCode (*reset)(MFN); | ||
33 | PetscErrorCode (*view)(MFN,PetscViewer); | ||
34 | }; | ||
35 | |||
36 | /* | ||
37 | Maximum number of monitors you can run with a single MFN | ||
38 | */ | ||
39 | #define MAXMFNMONITORS 5 | ||
40 | |||
41 | /* | ||
42 | Defines the MFN data structure. | ||
43 | */ | ||
44 | struct _p_MFN { | ||
45 | PETSCHEADER(struct _MFNOps); | ||
46 | /*------------------------- User parameters ---------------------------*/ | ||
47 | Mat A; /* the problem matrix */ | ||
48 | FN fn; /* which function to compute */ | ||
49 | PetscInt max_it; /* maximum number of iterations */ | ||
50 | PetscInt ncv; /* number of basis vectors */ | ||
51 | PetscReal tol; /* tolerance */ | ||
52 | PetscBool errorifnotconverged; /* error out if MFNSolve() does not converge */ | ||
53 | |||
54 | /*-------------- User-provided functions and contexts -----------------*/ | ||
55 | MFNMonitorFn *monitor[MAXMFNMONITORS]; | ||
56 | PetscCtxDestroyFn *monitordestroy[MAXMFNMONITORS]; | ||
57 | void *monitorcontext[MAXMFNMONITORS]; | ||
58 | PetscInt numbermonitors; | ||
59 | |||
60 | /*----------------- Child objects and working data -------------------*/ | ||
61 | BV V; /* set of basis vectors */ | ||
62 | Mat AT; /* the transpose of A, used in MFNSolveTranspose */ | ||
63 | PetscInt nwork; /* number of work vectors */ | ||
64 | Vec *work; /* work vectors */ | ||
65 | void *data; /* placeholder for solver-specific stuff */ | ||
66 | |||
67 | /* ----------------------- Status variables -------------------------- */ | ||
68 | PetscInt its; /* number of iterations so far computed */ | ||
69 | PetscInt nv; /* size of current Schur decomposition */ | ||
70 | PetscReal errest; /* error estimate */ | ||
71 | PetscReal bnorm; /* computed norm of right-hand side in current solve */ | ||
72 | PetscBool transpose_solve; /* solve transpose system instead */ | ||
73 | PetscInt setupcalled; | ||
74 | MFNConvergedReason reason; | ||
75 | }; | ||
76 | |||
77 | /* | ||
78 | MFN_CreateDenseMat - Creates a dense Mat of size k unless it already has that size | ||
79 | */ | ||
80 | 2486 | static inline PetscErrorCode MFN_CreateDenseMat(PetscInt k,Mat *A) | |
81 | { | ||
82 | 2486 | PetscBool create=PETSC_FALSE; | |
83 | 2486 | PetscInt m,n; | |
84 | |||
85 |
1/2✓ Branch 0 taken 4 times.
✗ Branch 1 not taken.
|
2486 | PetscFunctionBegin; |
86 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 20 times.
|
2486 | if (!*A) create=PETSC_TRUE; |
87 | else { | ||
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.
|
20 | PetscCall(MatGetSize(*A,&m,&n)); |
89 |
2/4✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 10 times.
|
20 | if (m!=k || n!=k) { |
90 | ✗ | PetscCall(MatDestroy(A)); | |
91 | create=PETSC_TRUE; | ||
92 | } | ||
93 | } | ||
94 |
4/6✓ Branch 0 taken 4 times.
✓ Branch 1 taken 16 times.
✓ Branch 2 taken 4 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 4 times.
|
2466 | if (create) PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,A)); |
95 |
6/12✓ Branch 0 taken 4 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 4 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 4 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 4 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 4 times.
|
485 | PetscFunctionReturn(PETSC_SUCCESS); |
96 | } | ||
97 | |||
98 | /* | ||
99 | MFN_CreateVec - Creates a Vec of size k unless it already has that size | ||
100 | */ | ||
101 | 2298 | static inline PetscErrorCode MFN_CreateVec(PetscInt k,Vec *v) | |
102 | { | ||
103 | 2298 | PetscBool create=PETSC_FALSE; | |
104 | 2298 | PetscInt n; | |
105 | |||
106 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
2298 | PetscFunctionBegin; |
107 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2298 | if (!*v) create=PETSC_TRUE; |
108 | else { | ||
109 |
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.
|
1334 | PetscCall(VecGetSize(*v,&n)); |
110 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
1334 | if (n!=k) { |
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.
|
1334 | PetscCall(VecDestroy(v)); |
112 | create=PETSC_TRUE; | ||
113 | } | ||
114 | } | ||
115 |
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.
|
2298 | if (create) PetscCall(VecCreateSeq(PETSC_COMM_SELF,k,v)); |
116 |
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.
|
449 | PetscFunctionReturn(PETSC_SUCCESS); |
117 | } | ||
118 |