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 | The ST interface routines, callable by users | ||
12 | */ | ||
13 | |||
14 | #include <slepc/private/stimpl.h> /*I "slepcst.h" I*/ | ||
15 | |||
16 | PetscClassId ST_CLASSID = 0; | ||
17 | PetscLogEvent ST_SetUp = 0,ST_ComputeOperator = 0,ST_Apply = 0,ST_ApplyTranspose = 0,ST_ApplyHermitianTranspose = 0,ST_MatSetUp = 0,ST_MatMult = 0,ST_MatMultTranspose = 0,ST_MatSolve = 0,ST_MatSolveTranspose = 0; | ||
18 | static PetscBool STPackageInitialized = PETSC_FALSE; | ||
19 | |||
20 | const char *STMatModes[] = {"COPY","INPLACE","SHELL","STMatMode","ST_MATMODE_",NULL}; | ||
21 | |||
22 | /*@C | ||
23 | STFinalizePackage - This function destroys everything in the Slepc interface | ||
24 | to the ST package. It is called from SlepcFinalize(). | ||
25 | |||
26 | Level: developer | ||
27 | |||
28 | .seealso: SlepcFinalize() | ||
29 | @*/ | ||
30 | 7999 | PetscErrorCode STFinalizePackage(void) | |
31 | { | ||
32 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
7999 | PetscFunctionBegin; |
33 |
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.
|
7999 | PetscCall(PetscFunctionListDestroy(&STList)); |
34 | 7999 | STPackageInitialized = PETSC_FALSE; | |
35 | 7999 | STRegisterAllCalled = PETSC_FALSE; | |
36 |
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.
|
7999 | PetscFunctionReturn(PETSC_SUCCESS); |
37 | } | ||
38 | |||
39 | /*@C | ||
40 | STInitializePackage - This function initializes everything in the ST package. | ||
41 | It is called from PetscDLLibraryRegister() when using dynamic libraries, and | ||
42 | on the first call to STCreate() when using static libraries. | ||
43 | |||
44 | Level: developer | ||
45 | |||
46 | .seealso: SlepcInitialize() | ||
47 | @*/ | ||
48 | 57187 | PetscErrorCode STInitializePackage(void) | |
49 | { | ||
50 | 57187 | char logList[256]; | |
51 | 57187 | PetscBool opt,pkg; | |
52 | 57187 | PetscClassId classids[1]; | |
53 | |||
54 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
57187 | PetscFunctionBegin; |
55 |
8/14✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 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.
|
57187 | if (STPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS); |
56 | 7999 | STPackageInitialized = PETSC_TRUE; | |
57 | /* Register Classes */ | ||
58 |
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.
|
7999 | PetscCall(PetscClassIdRegister("Spectral Transform",&ST_CLASSID)); |
59 | /* Register Constructors */ | ||
60 |
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.
|
7999 | PetscCall(STRegisterAll()); |
61 | /* Register Events */ | ||
62 |
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.
|
7999 | PetscCall(PetscLogEventRegister("STSetUp",ST_CLASSID,&ST_SetUp)); |
63 |
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.
|
7999 | PetscCall(PetscLogEventRegister("STComputeOperatr",ST_CLASSID,&ST_ComputeOperator)); |
64 |
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.
|
7999 | PetscCall(PetscLogEventRegister("STApply",ST_CLASSID,&ST_Apply)); |
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.
|
7999 | PetscCall(PetscLogEventRegister("STApplyTranspose",ST_CLASSID,&ST_ApplyTranspose)); |
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.
|
7999 | PetscCall(PetscLogEventRegister("STApplyHermTrans",ST_CLASSID,&ST_ApplyHermitianTranspose)); |
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.
|
7999 | PetscCall(PetscLogEventRegister("STMatSetUp",ST_CLASSID,&ST_MatSetUp)); |
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.
|
7999 | PetscCall(PetscLogEventRegister("STMatMult",ST_CLASSID,&ST_MatMult)); |
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.
|
7999 | PetscCall(PetscLogEventRegister("STMatMultTranspose",ST_CLASSID,&ST_MatMultTranspose)); |
70 |
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.
|
7999 | PetscCall(PetscLogEventRegister("STMatSolve",ST_CLASSID,&ST_MatSolve)); |
71 |
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.
|
7999 | PetscCall(PetscLogEventRegister("STMatSolveTranspose",ST_CLASSID,&ST_MatSolveTranspose)); |
72 | /* Process Info */ | ||
73 | 7999 | classids[0] = ST_CLASSID; | |
74 |
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.
|
7999 | PetscCall(PetscInfoProcessClass("st",1,&classids[0])); |
75 | /* Process summary exclusions */ | ||
76 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
7999 | PetscCall(PetscOptionsGetString(NULL,NULL,"-log_exclude",logList,sizeof(logList),&opt)); |
77 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
7999 | if (opt) { |
78 |
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.
|
84 | PetscCall(PetscStrInList("st",logList,',',&pkg)); |
79 |
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.
|
84 | if (pkg) PetscCall(PetscLogEventDeactivateClass(ST_CLASSID)); |
80 | } | ||
81 | /* Register package finalizer */ | ||
82 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
7999 | PetscCall(PetscRegisterFinalize(STFinalizePackage)); |
83 |
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.
|
1571 | PetscFunctionReturn(PETSC_SUCCESS); |
84 | } | ||
85 | |||
86 | /*@ | ||
87 | STReset - Resets the ST context to the initial state (prior to setup) | ||
88 | and destroys any allocated Vecs and Mats. | ||
89 | |||
90 | Collective | ||
91 | |||
92 | Input Parameter: | ||
93 | . st - the spectral transformation context | ||
94 | |||
95 | Level: advanced | ||
96 | |||
97 | .seealso: STDestroy() | ||
98 | @*/ | ||
99 | 20368 | PetscErrorCode STReset(ST st) | |
100 | { | ||
101 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
20368 | PetscFunctionBegin; |
102 |
3/14✓ 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.
✗ 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.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
|
20368 | if (st) PetscValidHeaderSpecific(st,ST_CLASSID,1); |
103 |
2/14✓ Branch 0 taken 8 times.
✓ Branch 1 taken 1 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.
|
18357 | if (!st) PetscFunctionReturn(PETSC_SUCCESS); |
104 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
20368 | STCheckNotSeized(st,1); |
105 |
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.
|
20368 | PetscTryTypeMethod(st,reset); |
106 |
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.
|
20368 | if (st->ksp) PetscCall(KSPReset(st->ksp)); |
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.
|
20368 | PetscCall(MatDestroyMatrices(PetscMax(2,st->nmat),&st->T)); |
108 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20368 | PetscCall(MatDestroyMatrices(PetscMax(2,st->nmat),&st->A)); |
109 | 20368 | st->nmat = 0; | |
110 |
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.
|
20368 | PetscCall(PetscFree(st->Astate)); |
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.
|
20368 | PetscCall(MatDestroy(&st->Op)); |
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.
|
20368 | PetscCall(MatDestroy(&st->P)); |
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.
|
20368 | PetscCall(MatDestroy(&st->Pmat)); |
114 |
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.
|
20368 | PetscCall(MatDestroyMatrices(st->nsplit,&st->Psplit)); |
115 | 20368 | st->nsplit = 0; | |
116 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20368 | PetscCall(VecDestroyVecs(st->nwork,&st->work)); |
117 | 20368 | st->nwork = 0; | |
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.
|
20368 | PetscCall(VecDestroy(&st->wb)); |
119 |
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.
|
20368 | PetscCall(VecDestroy(&st->wht)); |
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.
|
20368 | PetscCall(VecDestroy(&st->D)); |
121 | 20368 | st->state = ST_STATE_INITIAL; | |
122 | 20368 | st->opready = PETSC_FALSE; | |
123 |
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.
|
20368 | PetscFunctionReturn(PETSC_SUCCESS); |
124 | } | ||
125 | |||
126 | /*@ | ||
127 | STDestroy - Destroys ST context that was created with STCreate(). | ||
128 | |||
129 | Collective | ||
130 | |||
131 | Input Parameter: | ||
132 | . st - the spectral transformation context | ||
133 | |||
134 | Level: beginner | ||
135 | |||
136 | .seealso: STCreate(), STSetUp() | ||
137 | @*/ | ||
138 | 9393 | PetscErrorCode STDestroy(ST *st) | |
139 | { | ||
140 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
9393 | PetscFunctionBegin; |
141 |
8/14✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 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.
|
9393 | if (!*st) PetscFunctionReturn(PETSC_SUCCESS); |
142 |
2/12✗ 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 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
|
9189 | PetscValidHeaderSpecific(*st,ST_CLASSID,1); |
143 |
8/14✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 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.
|
9189 | if (--((PetscObject)*st)->refct > 0) { *st = NULL; PetscFunctionReturn(PETSC_SUCCESS); } |
144 |
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.
|
9148 | PetscCall(STReset(*st)); |
145 |
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.
|
9148 | PetscTryTypeMethod(*st,destroy); |
146 |
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.
|
9148 | PetscCall(KSPDestroy(&(*st)->ksp)); |
147 |
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.
|
9148 | PetscCall(PetscHeaderDestroy(st)); |
148 |
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.
|
1765 | PetscFunctionReturn(PETSC_SUCCESS); |
149 | } | ||
150 | |||
151 | /*@ | ||
152 | STCreate - Creates a spectral transformation context. | ||
153 | |||
154 | Collective | ||
155 | |||
156 | Input Parameter: | ||
157 | . comm - MPI communicator | ||
158 | |||
159 | Output Parameter: | ||
160 | . newst - location to put the spectral transformation context | ||
161 | |||
162 | Level: beginner | ||
163 | |||
164 | .seealso: STSetUp(), STApply(), STDestroy(), ST | ||
165 | @*/ | ||
166 | 9148 | PetscErrorCode STCreate(MPI_Comm comm,ST *newst) | |
167 | { | ||
168 | 9148 | ST st; | |
169 | |||
170 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
9148 | PetscFunctionBegin; |
171 |
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.
|
9148 | PetscAssertPointer(newst,2); |
172 |
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.
|
9148 | PetscCall(STInitializePackage()); |
173 |
7/12✓ 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 10 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 8 times.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
9148 | PetscCall(SlepcHeaderCreate(st,ST_CLASSID,"ST","Spectral Transformation","ST",comm,STDestroy,STView)); |
174 | |||
175 | 9148 | st->A = NULL; | |
176 | 9148 | st->nmat = 0; | |
177 | 9148 | st->sigma = 0.0; | |
178 | 9148 | st->defsigma = 0.0; | |
179 | 9148 | st->matmode = ST_MATMODE_COPY; | |
180 | 9148 | st->str = UNKNOWN_NONZERO_PATTERN; | |
181 | 9148 | st->transform = PETSC_FALSE; | |
182 | 9148 | st->structured = PETSC_FALSE; | |
183 | 9148 | st->D = NULL; | |
184 | 9148 | st->Pmat = NULL; | |
185 | 9148 | st->Pmat_set = PETSC_FALSE; | |
186 | 9148 | st->Psplit = NULL; | |
187 | 9148 | st->nsplit = 0; | |
188 | 9148 | st->strp = UNKNOWN_NONZERO_PATTERN; | |
189 | |||
190 | 9148 | st->ksp = NULL; | |
191 | 9148 | st->usesksp = PETSC_FALSE; | |
192 | 9148 | st->nwork = 0; | |
193 | 9148 | st->work = NULL; | |
194 | 9148 | st->wb = NULL; | |
195 | 9148 | st->wht = NULL; | |
196 | 9148 | st->state = ST_STATE_INITIAL; | |
197 | 9148 | st->Astate = NULL; | |
198 | 9148 | st->T = NULL; | |
199 | 9148 | st->Op = NULL; | |
200 | 9148 | st->opseized = PETSC_FALSE; | |
201 | 9148 | st->opready = PETSC_FALSE; | |
202 | 9148 | st->P = NULL; | |
203 | 9148 | st->M = NULL; | |
204 | 9148 | st->sigma_set = PETSC_FALSE; | |
205 | 9148 | st->asymm = PETSC_FALSE; | |
206 | 9148 | st->aherm = PETSC_FALSE; | |
207 | 9148 | st->data = NULL; | |
208 | |||
209 | 9148 | *newst = st; | |
210 |
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.
|
9148 | PetscFunctionReturn(PETSC_SUCCESS); |
211 | } | ||
212 | |||
213 | /* | ||
214 | Checks whether the ST matrices are all symmetric or hermitian. | ||
215 | */ | ||
216 | 9838 | static inline PetscErrorCode STMatIsSymmetricKnown(ST st,PetscBool *symm,PetscBool *herm) | |
217 | { | ||
218 | 9838 | PetscInt i; | |
219 | 9838 | PetscBool sbaij=PETSC_FALSE,set,flg=PETSC_FALSE; | |
220 | |||
221 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
9838 | PetscFunctionBegin; |
222 | /* check if problem matrices are all sbaij */ | ||
223 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
9886 | for (i=0;i<st->nmat;i++) { |
224 |
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.
|
9862 | PetscCall(PetscObjectTypeCompareAny((PetscObject)st->A[i],&sbaij,MATSEQSBAIJ,MATMPISBAIJ,"")); |
225 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
9862 | if (!sbaij) break; |
226 | } | ||
227 | /* check if user has set the symmetric flag */ | ||
228 | 9838 | *symm = PETSC_TRUE; | |
229 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
10496 | for (i=0;i<st->nmat;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.
|
10168 | PetscCall(MatIsSymmetricKnown(st->A[i],&set,&flg)); |
231 |
4/4✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
|
10168 | if (!set || !flg) { *symm = PETSC_FALSE; break; } |
232 | } | ||
233 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
9838 | if (sbaij) *symm = PETSC_TRUE; |
234 | #if defined(PETSC_USE_COMPLEX) | ||
235 | /* check if user has set the hermitian flag */ | ||
236 | 4850 | *herm = PETSC_TRUE; | |
237 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
5182 | for (i=0;i<st->nmat;i++) { |
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.
|
4994 | PetscCall(MatIsHermitianKnown(st->A[i],&set,&flg)); |
239 |
4/4✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 5 times.
✓ Branch 3 taken 5 times.
|
4994 | if (!set || !flg) { *herm = PETSC_FALSE; break; } |
240 | } | ||
241 | #else | ||
242 | 4988 | *herm = *symm; | |
243 | #endif | ||
244 |
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.
|
5938 | PetscFunctionReturn(PETSC_SUCCESS); |
245 | } | ||
246 | |||
247 | /*@ | ||
248 | STSetMatrices - Sets the matrices associated with the eigenvalue problem. | ||
249 | |||
250 | Collective | ||
251 | |||
252 | Input Parameters: | ||
253 | + st - the spectral transformation context | ||
254 | . n - number of matrices in array A | ||
255 | - A - the array of matrices associated with the eigensystem | ||
256 | |||
257 | Notes: | ||
258 | It must be called before STSetUp(). If it is called again after STSetUp() then | ||
259 | the ST object is reset. | ||
260 | |||
261 | Level: intermediate | ||
262 | |||
263 | .seealso: STGetMatrix(), STGetNumMatrices(), STSetUp(), STReset() | ||
264 | @*/ | ||
265 | 10425 | PetscErrorCode STSetMatrices(ST st,PetscInt n,Mat A[]) | |
266 | { | ||
267 | 10425 | PetscInt i; | |
268 | 10425 | PetscBool same=PETSC_TRUE; | |
269 | |||
270 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
10425 | PetscFunctionBegin; |
271 |
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.
|
10425 | PetscValidHeaderSpecific(st,ST_CLASSID,1); |
272 |
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.
|
10425 | PetscValidLogicalCollectiveInt(st,n,2); |
273 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
10425 | PetscCheck(n>0,PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_OUTOFRANGE,"Must have one or more matrices, you have %" PetscInt_FMT,n); |
274 |
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.
|
10425 | PetscAssertPointer(A,3); |
275 |
13/32✓ 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 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.
✗ Branch 28 not taken.
✓ Branch 29 taken 2 times.
✗ Branch 30 not taken.
✗ Branch 31 not taken.
|
10425 | PetscCheckSameComm(st,1,*A,3); |
276 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
10425 | STCheckNotSeized(st,1); |
277 |
3/6✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 10 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
10425 | PetscCheck(!st->nsplit || st->nsplit==n,PetscObjectComm((PetscObject)st),PETSC_ERR_SUP,"The number of matrices must be the same as in STSetSplitPreconditioner()"); |
278 | |||
279 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
10425 | if (st->state) { |
280 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
946 | if (n!=st->nmat) same = PETSC_FALSE; |
281 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2046 | for (i=0;same&&i<n;i++) { |
282 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1100 | if (A[i]!=st->A[i]) same = PETSC_FALSE; |
283 | } | ||
284 |
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.
|
946 | if (!same) PetscCall(STReset(st)); |
285 | } else same = PETSC_FALSE; | ||
286 |
2/2✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
|
160 | if (!same) { |
287 |
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.
|
9838 | PetscCall(MatDestroyMatrices(PetscMax(2,st->nmat),&st->A)); |
288 |
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.
|
9838 | PetscCall(PetscCalloc1(PetscMax(2,n),&st->A)); |
289 |
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.
|
9838 | PetscCall(PetscFree(st->Astate)); |
290 |
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.
|
9838 | PetscCall(PetscMalloc1(PetscMax(2,n),&st->Astate)); |
291 | } | ||
292 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
28811 | for (i=0;i<n;i++) { |
293 |
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.
|
18386 | PetscValidHeaderSpecific(A[i],MAT_CLASSID,3); |
294 |
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.
|
18386 | PetscCall(PetscObjectReference((PetscObject)A[i])); |
295 |
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.
|
18386 | PetscCall(MatDestroy(&st->A[i])); |
296 | 18386 | st->A[i] = A[i]; | |
297 | 18386 | st->Astate[i] = ((PetscObject)A[i])->state; | |
298 | } | ||
299 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
10425 | if (n==1) { |
300 | 6028 | st->A[1] = NULL; | |
301 | 6028 | st->Astate[1] = 0; | |
302 | } | ||
303 | 10425 | st->nmat = n; | |
304 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
10425 | if (same) st->state = ST_STATE_UPDATED; |
305 | 9838 | else st->state = ST_STATE_INITIAL; | |
306 |
3/6✓ Branch 0 taken 5 times.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 5 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
10425 | PetscCheck(!same || !st->Psplit,PetscObjectComm((PetscObject)st),PETSC_ERR_SUP,"Support for changing the matrices while using a split preconditioner is not implemented yet"); |
307 | 10425 | st->opready = PETSC_FALSE; | |
308 |
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.
|
10425 | if (!same) PetscCall(STMatIsSymmetricKnown(st,&st->asymm,&st->aherm)); |
309 |
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.
|
2017 | PetscFunctionReturn(PETSC_SUCCESS); |
310 | } | ||
311 | |||
312 | /*@ | ||
313 | STGetMatrix - Gets the matrices associated with the original eigensystem. | ||
314 | |||
315 | Collective | ||
316 | |||
317 | Input Parameters: | ||
318 | + st - the spectral transformation context | ||
319 | - k - the index of the requested matrix (starting in 0) | ||
320 | |||
321 | Output Parameters: | ||
322 | . A - the requested matrix | ||
323 | |||
324 | Level: intermediate | ||
325 | |||
326 | .seealso: STSetMatrices(), STGetNumMatrices() | ||
327 | @*/ | ||
328 | 111558 | PetscErrorCode STGetMatrix(ST st,PetscInt k,Mat *A) | |
329 | { | ||
330 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
111558 | PetscFunctionBegin; |
331 |
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.
|
111558 | PetscValidHeaderSpecific(st,ST_CLASSID,1); |
332 |
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.
|
111558 | PetscValidLogicalCollectiveInt(st,k,2); |
333 |
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.
|
111558 | PetscAssertPointer(A,3); |
334 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
111558 | STCheckMatrices(st,1); |
335 |
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.
|
111558 | PetscCheck(k>=0 && k<st->nmat,PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_OUTOFRANGE,"k must be between 0 and %" PetscInt_FMT,st->nmat-1); |
336 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
111558 | PetscCheck(((PetscObject)st->A[k])->state==st->Astate[k],PetscObjectComm((PetscObject)st),PETSC_ERR_SUP,"Cannot retrieve original matrices (have been modified)"); |
337 | 111558 | *A = st->A[k]; | |
338 |
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.
|
111558 | PetscFunctionReturn(PETSC_SUCCESS); |
339 | } | ||
340 | |||
341 | /*@ | ||
342 | STGetMatrixTransformed - Gets the matrices associated with the transformed eigensystem. | ||
343 | |||
344 | Not Collective | ||
345 | |||
346 | Input Parameters: | ||
347 | + st - the spectral transformation context | ||
348 | - k - the index of the requested matrix (starting in 0) | ||
349 | |||
350 | Output Parameters: | ||
351 | . T - the requested matrix | ||
352 | |||
353 | Level: developer | ||
354 | |||
355 | .seealso: STGetMatrix(), STGetNumMatrices() | ||
356 | @*/ | ||
357 | 2379 | PetscErrorCode STGetMatrixTransformed(ST st,PetscInt k,Mat *T) | |
358 | { | ||
359 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
2379 | PetscFunctionBegin; |
360 |
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.
|
2379 | PetscValidHeaderSpecific(st,ST_CLASSID,1); |
361 |
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.
|
2379 | PetscValidLogicalCollectiveInt(st,k,2); |
362 |
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.
|
2379 | PetscAssertPointer(T,3); |
363 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
2379 | STCheckMatrices(st,1); |
364 |
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.
|
2379 | PetscCheck(k>=0 && k<st->nmat,PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_OUTOFRANGE,"k must be between 0 and %" PetscInt_FMT,st->nmat-1); |
365 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
2379 | PetscCheck(st->T,PetscObjectComm((PetscObject)st),PETSC_ERR_POINTER,"There are no transformed matrices"); |
366 | 2379 | *T = st->T[k]; | |
367 |
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.
|
2379 | PetscFunctionReturn(PETSC_SUCCESS); |
368 | } | ||
369 | |||
370 | /*@ | ||
371 | STGetNumMatrices - Returns the number of matrices stored in the ST. | ||
372 | |||
373 | Not Collective | ||
374 | |||
375 | Input Parameter: | ||
376 | . st - the spectral transformation context | ||
377 | |||
378 | Output Parameters: | ||
379 | . n - the number of matrices passed in STSetMatrices() | ||
380 | |||
381 | Level: intermediate | ||
382 | |||
383 | .seealso: STSetMatrices() | ||
384 | @*/ | ||
385 | 72608 | PetscErrorCode STGetNumMatrices(ST st,PetscInt *n) | |
386 | { | ||
387 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
72608 | PetscFunctionBegin; |
388 |
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.
|
72608 | PetscValidHeaderSpecific(st,ST_CLASSID,1); |
389 |
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.
|
72608 | PetscAssertPointer(n,2); |
390 | 72608 | *n = st->nmat; | |
391 |
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.
|
72608 | PetscFunctionReturn(PETSC_SUCCESS); |
392 | } | ||
393 | |||
394 | /*@ | ||
395 | STResetMatrixState - Resets the stored state of the matrices in the ST. | ||
396 | |||
397 | Logically Collective | ||
398 | |||
399 | Input Parameter: | ||
400 | . st - the spectral transformation context | ||
401 | |||
402 | Note: | ||
403 | This is useful in solvers where the user matrices are modified during | ||
404 | the computation, as in nonlinear inverse iteration. The effect is that | ||
405 | STGetMatrix() will retrieve the modified matrices as if they were | ||
406 | the matrices originally provided by the user. | ||
407 | |||
408 | Level: developer | ||
409 | |||
410 | .seealso: STGetMatrix(), EPSPowerSetNonlinear() | ||
411 | @*/ | ||
412 | 14435 | PetscErrorCode STResetMatrixState(ST st) | |
413 | { | ||
414 | 14435 | PetscInt i; | |
415 | |||
416 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
14435 | PetscFunctionBegin; |
417 |
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.
|
14435 | PetscValidHeaderSpecific(st,ST_CLASSID,1); |
418 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
43305 | for (i=0;i<st->nmat;i++) st->Astate[i] = ((PetscObject)st->A[i])->state; |
419 |
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.
|
14435 | PetscFunctionReturn(PETSC_SUCCESS); |
420 | } | ||
421 | |||
422 | /*@ | ||
423 | STSetPreconditionerMat - Sets the matrix to be used to build the preconditioner. | ||
424 | |||
425 | Collective | ||
426 | |||
427 | Input Parameters: | ||
428 | + st - the spectral transformation context | ||
429 | - mat - the matrix that will be used in constructing the preconditioner | ||
430 | |||
431 | Notes: | ||
432 | This matrix will be passed to the internal KSP object (via the last argument | ||
433 | of KSPSetOperators()) as the matrix to be used when constructing the preconditioner. | ||
434 | If no matrix is set or mat is set to NULL, A-sigma*B will be used | ||
435 | to build the preconditioner, being sigma the value set by STSetShift(). | ||
436 | |||
437 | More precisely, this is relevant for spectral transformations that represent | ||
438 | a rational matrix function, and use a KSP object for the denominator, called | ||
439 | K in the description of STGetOperator(). It includes also the STPRECOND case. | ||
440 | If the user has a good approximation to matrix K that can be used to build a | ||
441 | cheap preconditioner, it can be passed with this function. Note that it affects | ||
442 | only the Pmat argument of KSPSetOperators(), not the Amat argument. | ||
443 | |||
444 | If a preconditioner matrix is set, the default is to use an iterative KSP | ||
445 | rather than a direct method. | ||
446 | |||
447 | An alternative to pass an approximation of A-sigma*B with this function is | ||
448 | to provide approximations of A and B via STSetSplitPreconditioner(). The | ||
449 | difference is that when sigma changes the preconditioner is recomputed. | ||
450 | |||
451 | Use NULL to remove a previously set matrix. | ||
452 | |||
453 | Level: advanced | ||
454 | |||
455 | .seealso: STGetPreconditionerMat(), STSetShift(), STGetOperator(), STSetSplitPreconditioner() | ||
456 | @*/ | ||
457 | 210 | PetscErrorCode STSetPreconditionerMat(ST st,Mat mat) | |
458 | { | ||
459 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
210 | PetscFunctionBegin; |
460 |
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.
|
210 | PetscValidHeaderSpecific(st,ST_CLASSID,1); |
461 |
2/2✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
|
210 | if (mat) { |
462 |
2/12✗ 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 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
|
208 | PetscValidHeaderSpecific(mat,MAT_CLASSID,2); |
463 |
13/32✓ 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 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.
✗ Branch 28 not taken.
✓ Branch 29 taken 2 times.
✗ Branch 30 not taken.
✗ Branch 31 not taken.
|
208 | PetscCheckSameComm(st,1,mat,2); |
464 | } | ||
465 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
210 | STCheckNotSeized(st,1); |
466 |
3/6✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 10 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
210 | PetscCheck(!mat || !st->Psplit,PetscObjectComm((PetscObject)st),PETSC_ERR_SUP,"Cannot call both STSetPreconditionerMat and STSetSplitPreconditioner"); |
467 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 5 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.
|
205 | if (mat) PetscCall(PetscObjectReference((PetscObject)mat)); |
468 |
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(MatDestroy(&st->Pmat)); |
469 | 210 | st->Pmat = mat; | |
470 | 210 | st->Pmat_set = mat? PETSC_TRUE: PETSC_FALSE; | |
471 | 210 | st->state = ST_STATE_INITIAL; | |
472 | 210 | st->opready = PETSC_FALSE; | |
473 |
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.
|
210 | PetscFunctionReturn(PETSC_SUCCESS); |
474 | } | ||
475 | |||
476 | /*@ | ||
477 | STGetPreconditionerMat - Returns the matrix previously set by STSetPreconditionerMat(). | ||
478 | |||
479 | Not Collective | ||
480 | |||
481 | Input Parameter: | ||
482 | . st - the spectral transformation context | ||
483 | |||
484 | Output Parameter: | ||
485 | . mat - the matrix that will be used in constructing the preconditioner or | ||
486 | NULL if no matrix was set by STSetPreconditionerMat(). | ||
487 | |||
488 | Level: advanced | ||
489 | |||
490 | .seealso: STSetPreconditionerMat() | ||
491 | @*/ | ||
492 | 460 | PetscErrorCode STGetPreconditionerMat(ST st,Mat *mat) | |
493 | { | ||
494 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
460 | PetscFunctionBegin; |
495 |
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.
|
460 | PetscValidHeaderSpecific(st,ST_CLASSID,1); |
496 |
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.
|
460 | PetscAssertPointer(mat,2); |
497 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
460 | *mat = st->Pmat_set? st->Pmat: NULL; |
498 |
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.
|
460 | PetscFunctionReturn(PETSC_SUCCESS); |
499 | } | ||
500 | |||
501 | /*@ | ||
502 | STSetSplitPreconditioner - Sets the matrices from which to build the preconditioner | ||
503 | in split form. | ||
504 | |||
505 | Collective | ||
506 | |||
507 | Input Parameters: | ||
508 | + st - the spectral transformation context | ||
509 | . n - number of matrices | ||
510 | . Psplit - array of matrices | ||
511 | - strp - structure flag for Psplit matrices | ||
512 | |||
513 | Notes: | ||
514 | The number of matrices passed here must be the same as in STSetMatrices(). | ||
515 | |||
516 | For linear eigenproblems, the preconditioner matrix is computed as | ||
517 | Pmat(sigma) = A0-sigma*B0, where A0 and B0 are approximations of A and B | ||
518 | (the eigenproblem matrices) provided via the Psplit array in this function. | ||
519 | Compared to STSetPreconditionerMat(), this function allows setting a preconditioner | ||
520 | in a way that is independent of the shift sigma. Whenever the value of sigma | ||
521 | changes the preconditioner is recomputed. | ||
522 | |||
523 | Similarly, for polynomial eigenproblems the matrix for the preconditioner | ||
524 | is expressed as Pmat(sigma) = sum_i Psplit_i*phi_i(sigma), for i=1,...,n, where | ||
525 | the phi_i's are the polynomial basis functions. | ||
526 | |||
527 | The structure flag provides information about the relative nonzero pattern of the | ||
528 | Psplit_i matrices, in the same way as in STSetMatStructure(). | ||
529 | |||
530 | Use n=0 to reset a previously set split preconditioner. | ||
531 | |||
532 | Level: advanced | ||
533 | |||
534 | .seealso: STGetSplitPreconditionerTerm(), STGetSplitPreconditionerInfo(), STSetPreconditionerMat(), STSetMatrices(), STSetMatStructure() | ||
535 | @*/ | ||
536 | 215 | PetscErrorCode STSetSplitPreconditioner(ST st,PetscInt n,Mat Psplit[],MatStructure strp) | |
537 | { | ||
538 | 215 | PetscInt i,N=0,M,M0=0,mloc,nloc,mloc0=0; | |
539 | |||
540 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
215 | PetscFunctionBegin; |
541 |
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.
|
215 | PetscValidHeaderSpecific(st,ST_CLASSID,1); |
542 |
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.
|
215 | PetscValidLogicalCollectiveInt(st,n,2); |
543 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
215 | PetscCheck(n>=0,PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_OUTOFRANGE,"Negative value of n = %" PetscInt_FMT,n); |
544 |
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.
|
215 | PetscCheck(!n || !st->Pmat_set,PetscObjectComm((PetscObject)st),PETSC_ERR_SUP,"Cannot call both STSetPreconditionerMat and STSetSplitPreconditioner"); |
545 |
6/8✓ Branch 0 taken 9 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 5 times.
✓ Branch 4 taken 3 times.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
215 | PetscCheck(!n || !st->nmat || st->nmat==n,PetscObjectComm((PetscObject)st),PETSC_ERR_SUP,"The number of matrices must be the same as in STSetMatrices()"); |
546 |
3/10✓ 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.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
|
43 | if (n) PetscAssertPointer(Psplit,3); |
547 |
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.
|
43 | PetscValidLogicalCollectiveEnum(st,strp,4); |
548 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
43 | STCheckNotSeized(st,1); |
549 | |||
550 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
700 | for (i=0;i<n;i++) { |
551 |
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.
|
485 | PetscValidHeaderSpecific(Psplit[i],MAT_CLASSID,3); |
552 |
13/32✓ 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 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.
✗ Branch 28 not taken.
✓ Branch 29 taken 2 times.
✗ Branch 30 not taken.
✗ Branch 31 not taken.
|
485 | PetscCheckSameComm(st,1,Psplit[i],3); |
553 |
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.
|
485 | PetscCall(MatGetSize(Psplit[i],&M,&N)); |
554 |
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.
|
485 | PetscCall(MatGetLocalSize(Psplit[i],&mloc,&nloc)); |
555 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
485 | PetscCheck(M==N,PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_WRONG,"Psplit[%" PetscInt_FMT "] is a non-square matrix (%" PetscInt_FMT " rows, %" PetscInt_FMT " cols)",i,M,N); |
556 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
485 | PetscCheck(mloc==nloc,PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_WRONG,"Psplit[%" PetscInt_FMT "] does not have equal row and column local sizes (%" PetscInt_FMT ", %" PetscInt_FMT ")",i,mloc,nloc); |
557 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
485 | if (!i) { M0 = M; mloc0 = mloc; } |
558 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
485 | PetscCheck(M==M0,PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_INCOMP,"Dimensions of Psplit[%" PetscInt_FMT "] do not match with previous matrices (%" PetscInt_FMT ", %" PetscInt_FMT ")",i,M,M0); |
559 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
485 | PetscCheck(mloc==mloc0,PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_INCOMP,"Local dimensions of Psplit[%" PetscInt_FMT "] do not match with previous matrices (%" PetscInt_FMT ", %" PetscInt_FMT ")",i,mloc,mloc0); |
560 |
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.
|
485 | PetscCall(PetscObjectReference((PetscObject)Psplit[i])); |
561 | } | ||
562 | |||
563 |
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.
|
215 | if (st->Psplit) PetscCall(MatDestroyMatrices(st->nsplit,&st->Psplit)); |
564 | |||
565 | /* allocate space and copy matrices */ | ||
566 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
215 | if (n) { |
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.
|
215 | PetscCall(PetscMalloc1(n,&st->Psplit)); |
568 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
700 | for (i=0;i<n;i++) st->Psplit[i] = Psplit[i]; |
569 | } | ||
570 | 215 | st->nsplit = n; | |
571 | 215 | st->strp = strp; | |
572 | 215 | st->state = ST_STATE_INITIAL; | |
573 |
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.
|
215 | PetscFunctionReturn(PETSC_SUCCESS); |
574 | } | ||
575 | |||
576 | /*@ | ||
577 | STGetSplitPreconditionerTerm - Gets the matrices associated with | ||
578 | the split preconditioner. | ||
579 | |||
580 | Not Collective | ||
581 | |||
582 | Input Parameters: | ||
583 | + st - the spectral transformation context | ||
584 | - k - the index of the requested matrix (starting in 0) | ||
585 | |||
586 | Output Parameter: | ||
587 | . Psplit - the returned matrix | ||
588 | |||
589 | Level: advanced | ||
590 | |||
591 | .seealso: STSetSplitPreconditioner(), STGetSplitPreconditionerInfo() | ||
592 | @*/ | ||
593 | 145 | PetscErrorCode STGetSplitPreconditionerTerm(ST st,PetscInt k,Mat *Psplit) | |
594 | { | ||
595 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
145 | PetscFunctionBegin; |
596 |
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.
|
145 | PetscValidHeaderSpecific(st,ST_CLASSID,1); |
597 |
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.
|
145 | PetscValidLogicalCollectiveInt(st,k,2); |
598 |
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.
|
145 | PetscAssertPointer(Psplit,3); |
599 |
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.
|
145 | PetscCheck(k>=0 && k<st->nsplit,PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_OUTOFRANGE,"k must be between 0 and %" PetscInt_FMT,st->nsplit-1); |
600 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
145 | PetscCheck(st->Psplit,PetscObjectComm((PetscObject)st),PETSC_ERR_ORDER,"You have not called STSetSplitPreconditioner()"); |
601 | 145 | *Psplit = st->Psplit[k]; | |
602 |
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.
|
145 | PetscFunctionReturn(PETSC_SUCCESS); |
603 | } | ||
604 | |||
605 | /*@ | ||
606 | STGetSplitPreconditionerInfo - Returns the number of matrices of the split | ||
607 | preconditioner, as well as the structure flag. | ||
608 | |||
609 | Not Collective | ||
610 | |||
611 | Input Parameter: | ||
612 | . st - the spectral transformation context | ||
613 | |||
614 | Output Parameters: | ||
615 | + n - the number of matrices passed in STSetSplitPreconditioner() | ||
616 | - strp - the matrix structure flag passed in STSetSplitPreconditioner() | ||
617 | |||
618 | Level: advanced | ||
619 | |||
620 | .seealso: STSetSplitPreconditioner(), STGetSplitPreconditionerTerm() | ||
621 | @*/ | ||
622 | 1379 | PetscErrorCode STGetSplitPreconditionerInfo(ST st,PetscInt *n,MatStructure *strp) | |
623 | { | ||
624 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
1379 | PetscFunctionBegin; |
625 |
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.
|
1379 | PetscValidHeaderSpecific(st,ST_CLASSID,1); |
626 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 5 times.
|
1379 | if (n) *n = st->nsplit; |
627 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1379 | if (strp) *strp = st->strp; |
628 |
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.
|
1379 | PetscFunctionReturn(PETSC_SUCCESS); |
629 | } | ||
630 | |||
631 | /*@ | ||
632 | STSetShift - Sets the shift associated with the spectral transformation. | ||
633 | |||
634 | Collective | ||
635 | |||
636 | Input Parameters: | ||
637 | + st - the spectral transformation context | ||
638 | - shift - the value of the shift | ||
639 | |||
640 | Notes: | ||
641 | In some spectral transformations, changing the shift may have associated | ||
642 | a lot of work, for example recomputing a factorization. | ||
643 | |||
644 | This function is normally not directly called by users, since the shift is | ||
645 | indirectly set by EPSSetTarget(). | ||
646 | |||
647 | Level: intermediate | ||
648 | |||
649 | .seealso: EPSSetTarget(), STGetShift(), STSetDefaultShift() | ||
650 | @*/ | ||
651 | 5990 | PetscErrorCode STSetShift(ST st,PetscScalar shift) | |
652 | { | ||
653 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
5990 | PetscFunctionBegin; |
654 |
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.
|
5990 | PetscValidHeaderSpecific(st,ST_CLASSID,1); |
655 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
5990 | PetscValidType(st,1); |
656 |
30/68✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✓ 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 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✓ Branch 21 taken 2 times.
✓ Branch 22 taken 2 times.
✗ Branch 23 not taken.
✗ Branch 24 not taken.
✓ Branch 25 taken 2 times.
✓ Branch 26 taken 2 times.
✗ Branch 27 not taken.
✗ Branch 28 not taken.
✓ Branch 29 taken 2 times.
✓ Branch 30 taken 2 times.
✗ Branch 31 not taken.
✗ Branch 32 not taken.
✓ Branch 33 taken 2 times.
✓ Branch 34 taken 2 times.
✗ Branch 35 not taken.
✓ Branch 36 taken 2 times.
✗ Branch 37 not taken.
✗ Branch 38 not taken.
✓ Branch 39 taken 2 times.
✗ Branch 40 not taken.
✗ Branch 41 not taken.
✗ Branch 42 not taken.
✓ Branch 43 taken 2 times.
✗ Branch 44 not taken.
✗ Branch 45 not taken.
✗ Branch 46 not taken.
✓ Branch 47 taken 2 times.
✓ Branch 48 taken 2 times.
✗ Branch 49 not taken.
✗ Branch 50 not taken.
✓ Branch 51 taken 2 times.
✓ Branch 52 taken 2 times.
✗ Branch 53 not taken.
✓ Branch 54 taken 2 times.
✗ Branch 55 not taken.
✗ Branch 56 not taken.
✓ Branch 57 taken 2 times.
✗ Branch 58 not taken.
✗ Branch 59 not taken.
✓ Branch 60 taken 2 times.
✗ Branch 61 not taken.
✓ Branch 62 taken 2 times.
✗ Branch 63 not taken.
✗ Branch 64 not taken.
✓ Branch 65 taken 2 times.
✗ Branch 66 not taken.
✗ Branch 67 not taken.
|
5990 | PetscValidLogicalCollectiveScalar(st,shift,2); |
657 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
5990 | if (st->sigma != shift) { |
658 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
5791 | STCheckNotSeized(st,1); |
659 |
7/10✓ 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 taken 8 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
|
5791 | if (st->state==ST_STATE_SETUP) PetscTryTypeMethod(st,setshift,shift); |
660 | 5791 | st->sigma = shift; | |
661 | } | ||
662 | 5990 | st->sigma_set = PETSC_TRUE; | |
663 |
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.
|
5990 | PetscFunctionReturn(PETSC_SUCCESS); |
664 | } | ||
665 | |||
666 | /*@ | ||
667 | STGetShift - Gets the shift associated with the spectral transformation. | ||
668 | |||
669 | Not Collective | ||
670 | |||
671 | Input Parameter: | ||
672 | . st - the spectral transformation context | ||
673 | |||
674 | Output Parameter: | ||
675 | . shift - the value of the shift | ||
676 | |||
677 | Level: intermediate | ||
678 | |||
679 | .seealso: STSetShift() | ||
680 | @*/ | ||
681 | 24434 | PetscErrorCode STGetShift(ST st,PetscScalar* shift) | |
682 | { | ||
683 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
24434 | PetscFunctionBegin; |
684 |
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.
|
24434 | PetscValidHeaderSpecific(st,ST_CLASSID,1); |
685 |
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.
|
24434 | PetscAssertPointer(shift,2); |
686 | 24434 | *shift = st->sigma; | |
687 |
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.
|
24434 | PetscFunctionReturn(PETSC_SUCCESS); |
688 | } | ||
689 | |||
690 | /*@ | ||
691 | STSetDefaultShift - Sets the value of the shift that should be employed if | ||
692 | the user did not specify one. | ||
693 | |||
694 | Logically Collective | ||
695 | |||
696 | Input Parameters: | ||
697 | + st - the spectral transformation context | ||
698 | - defaultshift - the default value of the shift | ||
699 | |||
700 | Level: developer | ||
701 | |||
702 | .seealso: STSetShift() | ||
703 | @*/ | ||
704 | 3303 | PetscErrorCode STSetDefaultShift(ST st,PetscScalar defaultshift) | |
705 | { | ||
706 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
3303 | PetscFunctionBegin; |
707 |
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.
|
3303 | PetscValidHeaderSpecific(st,ST_CLASSID,1); |
708 |
30/68✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✓ 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 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✓ Branch 21 taken 2 times.
✓ Branch 22 taken 2 times.
✗ Branch 23 not taken.
✗ Branch 24 not taken.
✓ Branch 25 taken 2 times.
✓ Branch 26 taken 2 times.
✗ Branch 27 not taken.
✗ Branch 28 not taken.
✓ Branch 29 taken 2 times.
✓ Branch 30 taken 2 times.
✗ Branch 31 not taken.
✗ Branch 32 not taken.
✓ Branch 33 taken 2 times.
✓ Branch 34 taken 2 times.
✗ Branch 35 not taken.
✓ Branch 36 taken 2 times.
✗ Branch 37 not taken.
✗ Branch 38 not taken.
✓ Branch 39 taken 2 times.
✗ Branch 40 not taken.
✗ Branch 41 not taken.
✗ Branch 42 not taken.
✓ Branch 43 taken 2 times.
✗ Branch 44 not taken.
✗ Branch 45 not taken.
✗ Branch 46 not taken.
✓ Branch 47 taken 2 times.
✓ Branch 48 taken 2 times.
✗ Branch 49 not taken.
✗ Branch 50 not taken.
✓ Branch 51 taken 2 times.
✓ Branch 52 taken 2 times.
✗ Branch 53 not taken.
✓ Branch 54 taken 2 times.
✗ Branch 55 not taken.
✗ Branch 56 not taken.
✓ Branch 57 taken 2 times.
✗ Branch 58 not taken.
✗ Branch 59 not taken.
✓ Branch 60 taken 2 times.
✗ Branch 61 not taken.
✓ Branch 62 taken 2 times.
✗ Branch 63 not taken.
✗ Branch 64 not taken.
✓ Branch 65 taken 2 times.
✗ Branch 66 not taken.
✗ Branch 67 not taken.
|
3303 | PetscValidLogicalCollectiveScalar(st,defaultshift,2); |
709 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
3303 | if (st->defsigma != defaultshift) { |
710 | 2109 | st->defsigma = defaultshift; | |
711 | 2109 | st->state = ST_STATE_INITIAL; | |
712 | 2109 | st->opready = PETSC_FALSE; | |
713 | } | ||
714 |
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.
|
3303 | PetscFunctionReturn(PETSC_SUCCESS); |
715 | } | ||
716 | |||
717 | /*@ | ||
718 | STScaleShift - Multiply the shift with a given factor. | ||
719 | |||
720 | Logically Collective | ||
721 | |||
722 | Input Parameters: | ||
723 | + st - the spectral transformation context | ||
724 | - factor - the scaling factor | ||
725 | |||
726 | Note: | ||
727 | This function does not update the transformation matrices, as opposed to | ||
728 | STSetShift(). | ||
729 | |||
730 | Level: developer | ||
731 | |||
732 | .seealso: STSetShift() | ||
733 | @*/ | ||
734 | 2208 | PetscErrorCode STScaleShift(ST st,PetscScalar factor) | |
735 | { | ||
736 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
2208 | PetscFunctionBegin; |
737 |
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.
|
2208 | PetscValidHeaderSpecific(st,ST_CLASSID,1); |
738 |
30/68✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✓ 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 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✓ Branch 21 taken 2 times.
✓ Branch 22 taken 2 times.
✗ Branch 23 not taken.
✗ Branch 24 not taken.
✓ Branch 25 taken 2 times.
✓ Branch 26 taken 2 times.
✗ Branch 27 not taken.
✗ Branch 28 not taken.
✓ Branch 29 taken 2 times.
✓ Branch 30 taken 2 times.
✗ Branch 31 not taken.
✗ Branch 32 not taken.
✓ Branch 33 taken 2 times.
✓ Branch 34 taken 2 times.
✗ Branch 35 not taken.
✓ Branch 36 taken 2 times.
✗ Branch 37 not taken.
✗ Branch 38 not taken.
✓ Branch 39 taken 2 times.
✗ Branch 40 not taken.
✗ Branch 41 not taken.
✗ Branch 42 not taken.
✓ Branch 43 taken 2 times.
✗ Branch 44 not taken.
✗ Branch 45 not taken.
✗ Branch 46 not taken.
✓ Branch 47 taken 2 times.
✓ Branch 48 taken 2 times.
✗ Branch 49 not taken.
✗ Branch 50 not taken.
✓ Branch 51 taken 2 times.
✓ Branch 52 taken 2 times.
✗ Branch 53 not taken.
✓ Branch 54 taken 2 times.
✗ Branch 55 not taken.
✗ Branch 56 not taken.
✓ Branch 57 taken 2 times.
✗ Branch 58 not taken.
✗ Branch 59 not taken.
✓ Branch 60 taken 2 times.
✗ Branch 61 not taken.
✓ Branch 62 taken 2 times.
✗ Branch 63 not taken.
✗ Branch 64 not taken.
✓ Branch 65 taken 2 times.
✗ Branch 66 not taken.
✗ Branch 67 not taken.
|
2208 | PetscValidLogicalCollectiveScalar(st,factor,2); |
739 | 2208 | st->sigma *= factor; | |
740 |
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.
|
2208 | PetscFunctionReturn(PETSC_SUCCESS); |
741 | } | ||
742 | |||
743 | /*@ | ||
744 | STSetBalanceMatrix - Sets the diagonal matrix to be used for balancing. | ||
745 | |||
746 | Collective | ||
747 | |||
748 | Input Parameters: | ||
749 | + st - the spectral transformation context | ||
750 | - D - the diagonal matrix (represented as a vector) | ||
751 | |||
752 | Notes: | ||
753 | If this matrix is set, STApply will effectively apply D*OP*D^{-1}. Use NULL | ||
754 | to reset a previously passed D. | ||
755 | |||
756 | Balancing is usually set via EPSSetBalance, but the advanced user may use | ||
757 | this function to bypass the usual balancing methods. | ||
758 | |||
759 | Level: developer | ||
760 | |||
761 | .seealso: EPSSetBalance(), STApply(), STGetBalanceMatrix() | ||
762 | @*/ | ||
763 | 9872 | PetscErrorCode STSetBalanceMatrix(ST st,Vec D) | |
764 | { | ||
765 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
9872 | PetscFunctionBegin; |
766 |
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.
|
9872 | PetscValidHeaderSpecific(st,ST_CLASSID,1); |
767 |
8/14✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 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.
|
9872 | if (st->D == D) PetscFunctionReturn(PETSC_SUCCESS); |
768 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
130 | STCheckNotSeized(st,1); |
769 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
130 | if (D) { |
770 |
2/12✗ 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 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
|
130 | PetscValidHeaderSpecific(D,VEC_CLASSID,2); |
771 |
13/32✓ 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 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.
✗ Branch 28 not taken.
✓ Branch 29 taken 2 times.
✗ Branch 30 not taken.
✗ Branch 31 not taken.
|
130 | PetscCheckSameComm(st,1,D,2); |
772 |
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.
|
130 | PetscCall(PetscObjectReference((PetscObject)D)); |
773 | } | ||
774 |
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.
|
130 | PetscCall(VecDestroy(&st->D)); |
775 | 130 | st->D = D; | |
776 | 130 | st->state = ST_STATE_INITIAL; | |
777 | 130 | st->opready = PETSC_FALSE; | |
778 |
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.
|
130 | PetscFunctionReturn(PETSC_SUCCESS); |
779 | } | ||
780 | |||
781 | /*@ | ||
782 | STGetBalanceMatrix - Gets the balance matrix used by the spectral transformation. | ||
783 | |||
784 | Not Collective | ||
785 | |||
786 | Input Parameter: | ||
787 | . st - the spectral transformation context | ||
788 | |||
789 | Output Parameter: | ||
790 | . D - the diagonal matrix (represented as a vector) | ||
791 | |||
792 | Note: | ||
793 | If the matrix was not set, a null pointer will be returned. | ||
794 | |||
795 | Level: developer | ||
796 | |||
797 | .seealso: STSetBalanceMatrix() | ||
798 | @*/ | ||
799 | ✗ | PetscErrorCode STGetBalanceMatrix(ST st,Vec *D) | |
800 | { | ||
801 | ✗ | PetscFunctionBegin; | |
802 | ✗ | PetscValidHeaderSpecific(st,ST_CLASSID,1); | |
803 | ✗ | PetscAssertPointer(D,2); | |
804 | ✗ | *D = st->D; | |
805 | ✗ | PetscFunctionReturn(PETSC_SUCCESS); | |
806 | } | ||
807 | |||
808 | /*@ | ||
809 | STMatCreateVecs - Get vector(s) compatible with the ST matrices. | ||
810 | |||
811 | Collective | ||
812 | |||
813 | Input Parameter: | ||
814 | . st - the spectral transformation context | ||
815 | |||
816 | Output Parameters: | ||
817 | + right - (optional) vector that the matrix can be multiplied against | ||
818 | - left - (optional) vector that the matrix vector product can be stored in | ||
819 | |||
820 | Level: developer | ||
821 | |||
822 | .seealso: STMatCreateVecsEmpty() | ||
823 | @*/ | ||
824 | 4227 | PetscErrorCode STMatCreateVecs(ST st,Vec *right,Vec *left) | |
825 | { | ||
826 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
4227 | PetscFunctionBegin; |
827 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
4227 | STCheckMatrices(st,1); |
828 |
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.
|
4227 | PetscCall(MatCreateVecs(st->A[0],right,left)); |
829 |
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.
|
814 | PetscFunctionReturn(PETSC_SUCCESS); |
830 | } | ||
831 | |||
832 | /*@ | ||
833 | STMatCreateVecsEmpty - Get vector(s) compatible with the ST matrices, i.e. with the same | ||
834 | parallel layout, but without internal array. | ||
835 | |||
836 | Collective | ||
837 | |||
838 | Input Parameter: | ||
839 | . st - the spectral transformation context | ||
840 | |||
841 | Output Parameters: | ||
842 | + right - (optional) vector that the matrix can be multiplied against | ||
843 | - left - (optional) vector that the matrix vector product can be stored in | ||
844 | |||
845 | Level: developer | ||
846 | |||
847 | .seealso: STMatCreateVecs(), MatCreateVecsEmpty() | ||
848 | @*/ | ||
849 | 9276 | PetscErrorCode STMatCreateVecsEmpty(ST st,Vec *right,Vec *left) | |
850 | { | ||
851 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
9276 | PetscFunctionBegin; |
852 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
9276 | STCheckMatrices(st,1); |
853 |
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.
|
9276 | PetscCall(MatCreateVecsEmpty(st->A[0],right,left)); |
854 |
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.
|
1812 | PetscFunctionReturn(PETSC_SUCCESS); |
855 | } | ||
856 | |||
857 | /*@ | ||
858 | STMatGetSize - Returns the number of rows and columns of the ST matrices. | ||
859 | |||
860 | Not Collective | ||
861 | |||
862 | Input Parameter: | ||
863 | . st - the spectral transformation context | ||
864 | |||
865 | Output Parameters: | ||
866 | + m - the number of global rows | ||
867 | - n - the number of global columns | ||
868 | |||
869 | Level: developer | ||
870 | |||
871 | .seealso: STMatGetLocalSize() | ||
872 | @*/ | ||
873 | 9742 | PetscErrorCode STMatGetSize(ST st,PetscInt *m,PetscInt *n) | |
874 | { | ||
875 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
9742 | PetscFunctionBegin; |
876 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
9742 | STCheckMatrices(st,1); |
877 |
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.
|
9742 | PetscCall(MatGetSize(st->A[0],m,n)); |
878 |
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.
|
1874 | PetscFunctionReturn(PETSC_SUCCESS); |
879 | } | ||
880 | |||
881 | /*@ | ||
882 | STMatGetLocalSize - Returns the number of local rows and columns of the ST matrices. | ||
883 | |||
884 | Not Collective | ||
885 | |||
886 | Input Parameter: | ||
887 | . st - the spectral transformation context | ||
888 | |||
889 | Output Parameters: | ||
890 | + m - the number of local rows | ||
891 | - n - the number of local columns | ||
892 | |||
893 | Level: developer | ||
894 | |||
895 | .seealso: STMatGetSize() | ||
896 | @*/ | ||
897 | 9742 | PetscErrorCode STMatGetLocalSize(ST st,PetscInt *m,PetscInt *n) | |
898 | { | ||
899 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
9742 | PetscFunctionBegin; |
900 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
9742 | STCheckMatrices(st,1); |
901 |
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.
|
9742 | PetscCall(MatGetLocalSize(st->A[0],m,n)); |
902 |
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.
|
1874 | PetscFunctionReturn(PETSC_SUCCESS); |
903 | } | ||
904 | |||
905 | /*@ | ||
906 | STSetOptionsPrefix - Sets the prefix used for searching for all | ||
907 | ST options in the database. | ||
908 | |||
909 | Logically Collective | ||
910 | |||
911 | Input Parameters: | ||
912 | + st - the spectral transformation context | ||
913 | - prefix - the prefix string to prepend to all ST option requests | ||
914 | |||
915 | Notes: | ||
916 | A hyphen (-) must NOT be given at the beginning of the prefix name. | ||
917 | The first character of all runtime options is AUTOMATICALLY the | ||
918 | hyphen. | ||
919 | |||
920 | Level: advanced | ||
921 | |||
922 | .seealso: STAppendOptionsPrefix(), STGetOptionsPrefix() | ||
923 | @*/ | ||
924 | 2413 | PetscErrorCode STSetOptionsPrefix(ST st,const char *prefix) | |
925 | { | ||
926 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
2413 | PetscFunctionBegin; |
927 |
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.
|
2413 | PetscValidHeaderSpecific(st,ST_CLASSID,1); |
928 |
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.
|
2413 | if (!st->ksp) PetscCall(STGetKSP(st,&st->ksp)); |
929 |
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.
|
2413 | PetscCall(KSPSetOptionsPrefix(st->ksp,prefix)); |
930 |
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.
|
2413 | PetscCall(KSPAppendOptionsPrefix(st->ksp,"st_")); |
931 |
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.
|
2413 | PetscCall(PetscObjectSetOptionsPrefix((PetscObject)st,prefix)); |
932 |
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.
|
464 | PetscFunctionReturn(PETSC_SUCCESS); |
933 | } | ||
934 | |||
935 | /*@ | ||
936 | STAppendOptionsPrefix - Appends to the prefix used for searching for all | ||
937 | ST options in the database. | ||
938 | |||
939 | Logically Collective | ||
940 | |||
941 | Input Parameters: | ||
942 | + st - the spectral transformation context | ||
943 | - prefix - the prefix string to prepend to all ST option requests | ||
944 | |||
945 | Notes: | ||
946 | A hyphen (-) must NOT be given at the beginning of the prefix name. | ||
947 | The first character of all runtime options is AUTOMATICALLY the | ||
948 | hyphen. | ||
949 | |||
950 | Level: advanced | ||
951 | |||
952 | .seealso: STSetOptionsPrefix(), STGetOptionsPrefix() | ||
953 | @*/ | ||
954 | 1934 | PetscErrorCode STAppendOptionsPrefix(ST st,const char *prefix) | |
955 | { | ||
956 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
1934 | PetscFunctionBegin; |
957 |
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.
|
1934 | PetscValidHeaderSpecific(st,ST_CLASSID,1); |
958 |
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.
|
1934 | PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)st,prefix)); |
959 |
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.
|
1934 | if (!st->ksp) PetscCall(STGetKSP(st,&st->ksp)); |
960 |
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.
|
1934 | PetscCall(KSPSetOptionsPrefix(st->ksp,((PetscObject)st)->prefix)); |
961 |
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.
|
1934 | PetscCall(KSPAppendOptionsPrefix(st->ksp,"st_")); |
962 |
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.
|
382 | PetscFunctionReturn(PETSC_SUCCESS); |
963 | } | ||
964 | |||
965 | /*@ | ||
966 | STGetOptionsPrefix - Gets the prefix used for searching for all | ||
967 | ST options in the database. | ||
968 | |||
969 | Not Collective | ||
970 | |||
971 | Input Parameters: | ||
972 | . st - the spectral transformation context | ||
973 | |||
974 | Output Parameters: | ||
975 | . prefix - pointer to the prefix string used, is returned | ||
976 | |||
977 | Note: | ||
978 | On the Fortran side, the user should pass in a string 'prefix' of | ||
979 | sufficient length to hold the prefix. | ||
980 | |||
981 | Level: advanced | ||
982 | |||
983 | .seealso: STSetOptionsPrefix(), STAppendOptionsPrefix() | ||
984 | @*/ | ||
985 | 21 | PetscErrorCode STGetOptionsPrefix(ST st,const char *prefix[]) | |
986 | { | ||
987 |
0/2✗ Branch 0 not taken.
✗ Branch 1 not taken.
|
21 | PetscFunctionBegin; |
988 |
0/16✗ Branch 0 not taken.
✗ Branch 1 not taken.
✗ 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.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
21 | PetscValidHeaderSpecific(st,ST_CLASSID,1); |
989 |
0/8✗ Branch 0 not taken.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
21 | PetscAssertPointer(prefix,2); |
990 |
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.
|
21 | PetscCall(PetscObjectGetOptionsPrefix((PetscObject)st,prefix)); |
991 | ✗ | PetscFunctionReturn(PETSC_SUCCESS); | |
992 | } | ||
993 | |||
994 | /*@ | ||
995 | STView - Prints the ST data structure. | ||
996 | |||
997 | Collective | ||
998 | |||
999 | Input Parameters: | ||
1000 | + st - the ST context | ||
1001 | - viewer - optional visualization context | ||
1002 | |||
1003 | Note: | ||
1004 | The available visualization contexts include | ||
1005 | + PETSC_VIEWER_STDOUT_SELF - standard output (default) | ||
1006 | - PETSC_VIEWER_STDOUT_WORLD - synchronized standard | ||
1007 | output where only the first processor opens | ||
1008 | the file. All other processors send their | ||
1009 | data to the first processor to print. | ||
1010 | |||
1011 | The user can open an alternative visualization contexts with | ||
1012 | PetscViewerASCIIOpen() (output to a specified file). | ||
1013 | |||
1014 | Level: beginner | ||
1015 | |||
1016 | .seealso: EPSView() | ||
1017 | @*/ | ||
1018 | 126 | PetscErrorCode STView(ST st,PetscViewer viewer) | |
1019 | { | ||
1020 | 126 | STType cstr; | |
1021 | 126 | char str[50]; | |
1022 | 126 | PetscBool isascii,isstring; | |
1023 | |||
1024 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
126 | PetscFunctionBegin; |
1025 |
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.
|
126 | PetscValidHeaderSpecific(st,ST_CLASSID,1); |
1026 |
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.
|
126 | if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)st),&viewer)); |
1027 |
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.
|
126 | PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); |
1028 |
13/32✓ 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 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.
✗ Branch 28 not taken.
✓ Branch 29 taken 2 times.
✗ Branch 30 not taken.
✗ Branch 31 not taken.
|
126 | PetscCheckSameComm(st,1,viewer,2); |
1029 | |||
1030 |
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.
|
126 | PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii)); |
1031 |
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.
|
126 | PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring)); |
1032 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
126 | if (isascii) { |
1033 |
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.
|
126 | PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)st,viewer)); |
1034 |
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.
|
126 | PetscCall(PetscViewerASCIIPushTab(viewer)); |
1035 |
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.
|
126 | PetscTryTypeMethod(st,view,viewer); |
1036 |
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.
|
126 | PetscCall(PetscViewerASCIIPopTab(viewer)); |
1037 |
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.
|
126 | PetscCall(SlepcSNPrintfScalar(str,sizeof(str),st->sigma,PETSC_FALSE)); |
1038 |
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.
|
126 | PetscCall(PetscViewerASCIIPrintf(viewer," shift: %s\n",str)); |
1039 |
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.
|
126 | PetscCall(PetscViewerASCIIPrintf(viewer," number of matrices: %" PetscInt_FMT "\n",st->nmat)); |
1040 |
2/3✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
|
126 | switch (st->matmode) { |
1041 | case ST_MATMODE_COPY: | ||
1042 | break; | ||
1043 | ✗ | case ST_MATMODE_INPLACE: | |
1044 | ✗ | PetscCall(PetscViewerASCIIPrintf(viewer," shifting the matrix and unshifting at exit\n")); | |
1045 | break; | ||
1046 | 10 | case ST_MATMODE_SHELL: | |
1047 |
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.
|
10 | PetscCall(PetscViewerASCIIPrintf(viewer," using a shell matrix\n")); |
1048 | break; | ||
1049 | } | ||
1050 |
7/10✓ 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 taken 8 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
|
126 | if (st->nmat>1 && st->matmode != ST_MATMODE_SHELL) PetscCall(PetscViewerASCIIPrintf(viewer," nonzero pattern of the matrices: %s\n",MatStructures[st->str])); |
1051 |
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.
|
126 | if (st->Psplit) PetscCall(PetscViewerASCIIPrintf(viewer," using split preconditioner matrices with %s\n",MatStructures[st->strp])); |
1052 |
3/10✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 10 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
|
126 | if (st->transform && st->nmat>2) PetscCall(PetscViewerASCIIPrintf(viewer," computing transformed matrices\n")); |
1053 |
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.
|
126 | if (st->structured) PetscCall(PetscViewerASCIIPrintf(viewer," exploiting structure in the application of the operator\n")); |
1054 | ✗ | } else if (isstring) { | |
1055 | ✗ | PetscCall(STGetType(st,&cstr)); | |
1056 | ✗ | PetscCall(PetscViewerStringSPrintf(viewer," %-7.7s",cstr)); | |
1057 | ✗ | PetscTryTypeMethod(st,view,viewer); | |
1058 | } | ||
1059 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
126 | if (st->usesksp) { |
1060 |
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.
|
59 | if (!st->ksp) PetscCall(STGetKSP(st,&st->ksp)); |
1061 |
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.
|
59 | PetscCall(PetscViewerASCIIPushTab(viewer)); |
1062 |
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.
|
59 | PetscCall(KSPView(st->ksp,viewer)); |
1063 |
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.
|
59 | PetscCall(PetscViewerASCIIPopTab(viewer)); |
1064 | } | ||
1065 |
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.
|
24 | PetscFunctionReturn(PETSC_SUCCESS); |
1066 | } | ||
1067 | |||
1068 | /*@ | ||
1069 | STViewFromOptions - View from options | ||
1070 | |||
1071 | Collective | ||
1072 | |||
1073 | Input Parameters: | ||
1074 | + st - the spectral transformation context | ||
1075 | . obj - optional object | ||
1076 | - name - command line option | ||
1077 | |||
1078 | Level: intermediate | ||
1079 | |||
1080 | .seealso: STView(), STCreate() | ||
1081 | @*/ | ||
1082 | ✗ | PetscErrorCode STViewFromOptions(ST st,PetscObject obj,const char name[]) | |
1083 | { | ||
1084 | ✗ | PetscFunctionBegin; | |
1085 | ✗ | PetscValidHeaderSpecific(st,ST_CLASSID,1); | |
1086 | ✗ | PetscCall(PetscObjectViewFromOptions((PetscObject)st,obj,name)); | |
1087 | ✗ | PetscFunctionReturn(PETSC_SUCCESS); | |
1088 | } | ||
1089 | |||
1090 | /*@C | ||
1091 | STRegister - Adds a method to the spectral transformation package. | ||
1092 | |||
1093 | Not Collective | ||
1094 | |||
1095 | Input Parameters: | ||
1096 | + name - name of a new user-defined transformation | ||
1097 | - function - routine to create method context | ||
1098 | |||
1099 | Notes: | ||
1100 | STRegister() may be called multiple times to add several user-defined | ||
1101 | spectral transformations. | ||
1102 | |||
1103 | Example Usage: | ||
1104 | .vb | ||
1105 | STRegister("my_transform",MyTransformCreate); | ||
1106 | .ve | ||
1107 | |||
1108 | Then, your spectral transform can be chosen with the procedural interface via | ||
1109 | $ STSetType(st,"my_transform") | ||
1110 | or at runtime via the option | ||
1111 | $ -st_type my_transform | ||
1112 | |||
1113 | Level: advanced | ||
1114 | |||
1115 | .seealso: STRegisterAll() | ||
1116 | @*/ | ||
1117 | 47994 | PetscErrorCode STRegister(const char *name,PetscErrorCode (*function)(ST)) | |
1118 | { | ||
1119 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
47994 | PetscFunctionBegin; |
1120 |
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.
|
47994 | PetscCall(STInitializePackage()); |
1121 |
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.
|
47994 | PetscCall(PetscFunctionListAdd(&STList,name,function)); |
1122 |
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.
|
9426 | PetscFunctionReturn(PETSC_SUCCESS); |
1123 | } | ||
1124 |