| 1887 |
jroman |
1 |
/*
|
|
|
2 |
QEP routines related to options that can be set via the command-line
|
|
|
3 |
or procedurally.
|
|
|
4 |
|
|
|
5 |
- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
|
|
|
6 |
SLEPc - Scalable Library for Eigenvalue Problem Computations
|
| 2116 |
eromero |
7 |
Copyright (c) 2002-2010, Universidad Politecnica de Valencia, Spain
|
| 1887 |
jroman |
8 |
|
|
|
9 |
This file is part of SLEPc.
|
|
|
10 |
|
|
|
11 |
SLEPc is free software: you can redistribute it and/or modify it under the
|
|
|
12 |
terms of version 3 of the GNU Lesser General Public License as published by
|
|
|
13 |
the Free Software Foundation.
|
|
|
14 |
|
|
|
15 |
SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
|
|
|
16 |
WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
|
|
|
17 |
FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
|
|
|
18 |
more details.
|
|
|
19 |
|
|
|
20 |
You should have received a copy of the GNU Lesser General Public License
|
|
|
21 |
along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
|
|
|
22 |
- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
|
|
|
23 |
*/
|
|
|
24 |
|
| 2284 |
jroman |
25 |
#include <private/qepimpl.h> /*I "slepcqep.h" I*/
|
| 1887 |
jroman |
26 |
|
|
|
27 |
#undef __FUNCT__
|
|
|
28 |
#define __FUNCT__ "QEPSetFromOptions"
|
|
|
29 |
/*@
|
|
|
30 |
QEPSetFromOptions - Sets QEP options from the options database.
|
|
|
31 |
This routine must be called before QEPSetUp() if the user is to be
|
|
|
32 |
allowed to set the solver type.
|
|
|
33 |
|
|
|
34 |
Collective on QEP
|
|
|
35 |
|
|
|
36 |
Input Parameters:
|
|
|
37 |
. qep - the quadratic eigensolver context
|
|
|
38 |
|
|
|
39 |
Notes:
|
|
|
40 |
To see all options, run your program with the -help option.
|
|
|
41 |
|
|
|
42 |
Level: beginner
|
|
|
43 |
@*/
|
|
|
44 |
PetscErrorCode QEPSetFromOptions(QEP qep)
|
|
|
45 |
{
|
| 2317 |
jroman |
46 |
PetscErrorCode ierr;
|
|
|
47 |
char type[256],monfilename[PETSC_MAX_PATH_LEN];
|
|
|
48 |
PetscBool flg,val;
|
|
|
49 |
PetscReal r;
|
|
|
50 |
PetscInt i,j,k;
|
| 1887 |
jroman |
51 |
PetscViewerASCIIMonitor monviewer;
|
| 2317 |
jroman |
52 |
SlepcConvMonitor ctx;
|
| 1887 |
jroman |
53 |
|
|
|
54 |
PetscFunctionBegin;
|
| 2213 |
jroman |
55 |
PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
|
| 1887 |
jroman |
56 |
ierr = PetscOptionsBegin(((PetscObject)qep)->comm,((PetscObject)qep)->prefix,"Quadratic Eigenvalue Problem (QEP) Solver Options","QEP");CHKERRQ(ierr);
|
|
|
57 |
ierr = PetscOptionsList("-qep_type","Quadratic Eigenvalue Problem method","QEPSetType",QEPList,(char*)(((PetscObject)qep)->type_name?((PetscObject)qep)->type_name:QEPLINEAR),type,256,&flg);CHKERRQ(ierr);
|
|
|
58 |
if (flg) {
|
|
|
59 |
ierr = QEPSetType(qep,type);CHKERRQ(ierr);
|
|
|
60 |
} else if (!((PetscObject)qep)->type_name) {
|
|
|
61 |
ierr = QEPSetType(qep,QEPLINEAR);CHKERRQ(ierr);
|
|
|
62 |
}
|
|
|
63 |
|
| 2216 |
jroman |
64 |
ierr = PetscOptionsBoolGroupBegin("-qep_general","general quadratic eigenvalue problem","QEPSetProblemType",&flg);CHKERRQ(ierr);
|
| 1907 |
jroman |
65 |
if (flg) {ierr = QEPSetProblemType(qep,QEP_GENERAL);CHKERRQ(ierr);}
|
| 2216 |
jroman |
66 |
ierr = PetscOptionsBoolGroup("-qep_hermitian","hermitian quadratic eigenvalue problem","QEPSetProblemType",&flg);CHKERRQ(ierr);
|
| 1907 |
jroman |
67 |
if (flg) {ierr = QEPSetProblemType(qep,QEP_HERMITIAN);CHKERRQ(ierr);}
|
| 2216 |
jroman |
68 |
ierr = PetscOptionsBoolGroupEnd("-qep_gyroscopic","gyroscopic quadratic eigenvalue problem","QEPSetProblemType",&flg);CHKERRQ(ierr);
|
| 1907 |
jroman |
69 |
if (flg) {ierr = QEPSetProblemType(qep,QEP_GYROSCOPIC);CHKERRQ(ierr);}
|
|
|
70 |
|
| 1918 |
jroman |
71 |
r = PETSC_IGNORE;
|
|
|
72 |
ierr = PetscOptionsReal("-qep_scale","Scale factor","QEPSetScaleFactor",qep->sfactor,&r,PETSC_NULL);CHKERRQ(ierr);
|
|
|
73 |
ierr = QEPSetScaleFactor(qep,r);CHKERRQ(ierr);
|
|
|
74 |
|
| 1887 |
jroman |
75 |
r = i = PETSC_IGNORE;
|
|
|
76 |
ierr = PetscOptionsInt("-qep_max_it","Maximum number of iterations","QEPSetTolerances",qep->max_it,&i,PETSC_NULL);CHKERRQ(ierr);
|
|
|
77 |
ierr = PetscOptionsReal("-qep_tol","Tolerance","QEPSetTolerances",qep->tol,&r,PETSC_NULL);CHKERRQ(ierr);
|
|
|
78 |
ierr = QEPSetTolerances(qep,r,i);CHKERRQ(ierr);
|
| 2216 |
jroman |
79 |
ierr = PetscOptionsBoolGroupBegin("-qep_convergence_default","Default (relative error) convergence test","QEPSetConvergenceTest",&flg);CHKERRQ(ierr);
|
| 1887 |
jroman |
80 |
if (flg) {ierr = QEPSetConvergenceTest(qep,QEPDefaultConverged,PETSC_NULL);CHKERRQ(ierr);}
|
| 2216 |
jroman |
81 |
ierr = PetscOptionsBoolGroupEnd("-qep_convergence_absolute","Absolute error convergence test","QEPSetConvergenceTest",&flg);CHKERRQ(ierr);
|
| 1887 |
jroman |
82 |
if (flg) {ierr = QEPSetConvergenceTest(qep,QEPAbsoluteConverged,PETSC_NULL);CHKERRQ(ierr);}
|
|
|
83 |
|
|
|
84 |
i = j = k = PETSC_IGNORE;
|
|
|
85 |
ierr = PetscOptionsInt("-qep_nev","Number of eigenvalues to compute","QEPSetDimensions",qep->nev,&i,PETSC_NULL);CHKERRQ(ierr);
|
|
|
86 |
ierr = PetscOptionsInt("-qep_ncv","Number of basis vectors","QEPSetDimensions",qep->ncv,&j,PETSC_NULL);CHKERRQ(ierr);
|
|
|
87 |
ierr = PetscOptionsInt("-qep_mpd","Maximum dimension of projected problem","QEPSetDimensions",qep->mpd,&k,PETSC_NULL);CHKERRQ(ierr);
|
|
|
88 |
ierr = QEPSetDimensions(qep,i,j,k);CHKERRQ(ierr);
|
|
|
89 |
|
|
|
90 |
/* -----------------------------------------------------------------------*/
|
|
|
91 |
/*
|
|
|
92 |
Cancels all monitors hardwired into code before call to QEPSetFromOptions()
|
|
|
93 |
*/
|
|
|
94 |
flg = PETSC_FALSE;
|
| 2216 |
jroman |
95 |
ierr = PetscOptionsBool("-qep_monitor_cancel","Remove any hardwired monitor routines","QEPMonitorCancel",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
|
| 1887 |
jroman |
96 |
if (flg) {
|
| 2330 |
jroman |
97 |
ierr = QEPMonitorCancel(qep);CHKERRQ(ierr);
|
| 1887 |
jroman |
98 |
}
|
|
|
99 |
/*
|
|
|
100 |
Prints approximate eigenvalues and error estimates at each iteration
|
|
|
101 |
*/
|
| 2054 |
eromero |
102 |
ierr = PetscOptionsString("-qep_monitor","Monitor first unconverged approximate eigenvalue and error estimate","QEPMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
|
| 1887 |
jroman |
103 |
if (flg) {
|
|
|
104 |
ierr = PetscViewerASCIIMonitorCreate(((PetscObject)qep)->comm,monfilename,((PetscObject)qep)->tablevel,&monviewer);CHKERRQ(ierr);
|
| 2054 |
eromero |
105 |
ierr = QEPMonitorSet(qep,QEPMonitorFirst,monviewer,(PetscErrorCode (*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
|
| 1887 |
jroman |
106 |
}
|
|
|
107 |
ierr = PetscOptionsString("-qep_monitor_conv","Monitor approximate eigenvalues and error estimates as they converge","QEPMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
|
|
|
108 |
if (flg) {
|
| 2311 |
jroman |
109 |
ierr = PetscNew(struct _n_SlepcConvMonitor,&ctx);CHKERRQ(ierr);
|
| 1951 |
jroman |
110 |
ierr = PetscViewerASCIIMonitorCreate(((PetscObject)qep)->comm,monfilename,((PetscObject)qep)->tablevel,&ctx->viewer);CHKERRQ(ierr);
|
| 2314 |
jroman |
111 |
ierr = QEPMonitorSet(qep,QEPMonitorConverged,ctx,(PetscErrorCode (*)(void*))SlepcConvMonitorDestroy);CHKERRQ(ierr);
|
| 1887 |
jroman |
112 |
}
|
| 2054 |
eromero |
113 |
ierr = PetscOptionsString("-qep_monitor_all","Monitor approximate eigenvalues and error estimates","QEPMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
|
| 1887 |
jroman |
114 |
if (flg) {
|
|
|
115 |
ierr = PetscViewerASCIIMonitorCreate(((PetscObject)qep)->comm,monfilename,((PetscObject)qep)->tablevel,&monviewer);CHKERRQ(ierr);
|
| 2054 |
eromero |
116 |
ierr = QEPMonitorSet(qep,QEPMonitorAll,monviewer,(PetscErrorCode (*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
|
|
|
117 |
ierr = QEPSetTrackAll(qep,PETSC_TRUE);CHKERRQ(ierr);
|
| 1887 |
jroman |
118 |
}
|
|
|
119 |
flg = PETSC_FALSE;
|
| 2216 |
jroman |
120 |
ierr = PetscOptionsBool("-qep_monitor_draw","Monitor first unconverged approximate error estimate graphically","QEPMonitorSet",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
|
| 1887 |
jroman |
121 |
if (flg) {
|
|
|
122 |
ierr = QEPMonitorSet(qep,QEPMonitorLG,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
|
|
|
123 |
}
|
| 2054 |
eromero |
124 |
flg = PETSC_FALSE;
|
| 2216 |
jroman |
125 |
ierr = PetscOptionsBool("-qep_monitor_draw_all","Monitor error estimates graphically","QEPMonitorSet",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
|
| 2054 |
eromero |
126 |
if (flg) {
|
|
|
127 |
ierr = QEPMonitorSet(qep,QEPMonitorLGAll,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
|
|
|
128 |
ierr = QEPSetTrackAll(qep,PETSC_TRUE);CHKERRQ(ierr);
|
|
|
129 |
}
|
| 1887 |
jroman |
130 |
/* -----------------------------------------------------------------------*/
|
|
|
131 |
|
| 2216 |
jroman |
132 |
ierr = PetscOptionsBoolGroupBegin("-qep_largest_magnitude","compute largest eigenvalues in magnitude","QEPSetWhichEigenpairs",&flg);CHKERRQ(ierr);
|
| 1887 |
jroman |
133 |
if (flg) {ierr = QEPSetWhichEigenpairs(qep,QEP_LARGEST_MAGNITUDE);CHKERRQ(ierr);}
|
| 2216 |
jroman |
134 |
ierr = PetscOptionsBoolGroup("-qep_smallest_magnitude","compute smallest eigenvalues in magnitude","QEPSetWhichEigenpairs",&flg);CHKERRQ(ierr);
|
| 1887 |
jroman |
135 |
if (flg) {ierr = QEPSetWhichEigenpairs(qep,QEP_SMALLEST_MAGNITUDE);CHKERRQ(ierr);}
|
| 2216 |
jroman |
136 |
ierr = PetscOptionsBoolGroup("-qep_largest_real","compute largest real parts","QEPSetWhichEigenpairs",&flg);CHKERRQ(ierr);
|
| 1887 |
jroman |
137 |
if (flg) {ierr = QEPSetWhichEigenpairs(qep,QEP_LARGEST_REAL);CHKERRQ(ierr);}
|
| 2216 |
jroman |
138 |
ierr = PetscOptionsBoolGroup("-qep_smallest_real","compute smallest real parts","QEPSetWhichEigenpairs",&flg);CHKERRQ(ierr);
|
| 1887 |
jroman |
139 |
if (flg) {ierr = QEPSetWhichEigenpairs(qep,QEP_SMALLEST_REAL);CHKERRQ(ierr);}
|
| 2216 |
jroman |
140 |
ierr = PetscOptionsBoolGroup("-qep_largest_imaginary","compute largest imaginary parts","QEPSetWhichEigenpairs",&flg);CHKERRQ(ierr);
|
| 1887 |
jroman |
141 |
if (flg) {ierr = QEPSetWhichEigenpairs(qep,QEP_LARGEST_IMAGINARY);CHKERRQ(ierr);}
|
| 2216 |
jroman |
142 |
ierr = PetscOptionsBoolGroupEnd("-qep_smallest_imaginary","compute smallest imaginary parts","QEPSetWhichEigenpairs",&flg);CHKERRQ(ierr);
|
| 1887 |
jroman |
143 |
if (flg) {ierr = QEPSetWhichEigenpairs(qep,QEP_SMALLEST_IMAGINARY);CHKERRQ(ierr);}
|
|
|
144 |
|
| 2216 |
jroman |
145 |
ierr = PetscOptionsBool("-qep_left_vectors","Compute left eigenvectors also","QEPSetLeftVectorsWanted",qep->leftvecs,&val,&flg);CHKERRQ(ierr);
|
| 1944 |
jroman |
146 |
if (flg) {
|
|
|
147 |
ierr = QEPSetLeftVectorsWanted(qep,val);CHKERRQ(ierr);
|
|
|
148 |
}
|
|
|
149 |
|
| 1887 |
jroman |
150 |
ierr = PetscOptionsName("-qep_view","Print detailed information on solver used","QEPView",0);CHKERRQ(ierr);
|
|
|
151 |
ierr = PetscOptionsName("-qep_view_binary","Save the matrices associated to the eigenproblem","QEPSetFromOptions",0);CHKERRQ(ierr);
|
|
|
152 |
ierr = PetscOptionsName("-qep_plot_eigs","Make a plot of the computed eigenvalues","QEPSolve",0);CHKERRQ(ierr);
|
|
|
153 |
|
|
|
154 |
if (qep->ops->setfromoptions) {
|
|
|
155 |
ierr = (*qep->ops->setfromoptions)(qep);CHKERRQ(ierr);
|
|
|
156 |
}
|
|
|
157 |
ierr = PetscOptionsEnd();CHKERRQ(ierr);
|
|
|
158 |
|
| 2330 |
jroman |
159 |
ierr = IPSetFromOptions(qep->ip);CHKERRQ(ierr);
|
| 1887 |
jroman |
160 |
PetscFunctionReturn(0);
|
|
|
161 |
}
|
|
|
162 |
|
|
|
163 |
#undef __FUNCT__
|
|
|
164 |
#define __FUNCT__ "QEPGetTolerances"
|
|
|
165 |
/*@
|
|
|
166 |
QEPGetTolerances - Gets the tolerance and maximum iteration count used
|
|
|
167 |
by the QEP convergence tests.
|
|
|
168 |
|
|
|
169 |
Not Collective
|
|
|
170 |
|
|
|
171 |
Input Parameter:
|
|
|
172 |
. qep - the quadratic eigensolver context
|
|
|
173 |
|
|
|
174 |
Output Parameters:
|
|
|
175 |
+ tol - the convergence tolerance
|
|
|
176 |
- maxits - maximum number of iterations
|
|
|
177 |
|
|
|
178 |
Notes:
|
|
|
179 |
The user can specify PETSC_NULL for any parameter that is not needed.
|
|
|
180 |
|
|
|
181 |
Level: intermediate
|
|
|
182 |
|
|
|
183 |
.seealso: QEPSetTolerances()
|
|
|
184 |
@*/
|
|
|
185 |
PetscErrorCode QEPGetTolerances(QEP qep,PetscReal *tol,PetscInt *maxits)
|
|
|
186 |
{
|
|
|
187 |
PetscFunctionBegin;
|
| 2213 |
jroman |
188 |
PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
|
| 1887 |
jroman |
189 |
if (tol) *tol = qep->tol;
|
|
|
190 |
if (maxits) *maxits = qep->max_it;
|
|
|
191 |
PetscFunctionReturn(0);
|
|
|
192 |
}
|
|
|
193 |
|
|
|
194 |
#undef __FUNCT__
|
| 1918 |
jroman |
195 |
#define __FUNCT__ "QEPSetTolerances"
|
| 1887 |
jroman |
196 |
/*@
|
|
|
197 |
QEPSetTolerances - Sets the tolerance and maximum iteration count used
|
|
|
198 |
by the QEP convergence tests.
|
|
|
199 |
|
| 2328 |
jroman |
200 |
Logically Collective on QEP
|
| 1887 |
jroman |
201 |
|
|
|
202 |
Input Parameters:
|
|
|
203 |
+ qep - the quadratic eigensolver context
|
|
|
204 |
. tol - the convergence tolerance
|
|
|
205 |
- maxits - maximum number of iterations to use
|
|
|
206 |
|
|
|
207 |
Options Database Keys:
|
|
|
208 |
+ -qep_tol <tol> - Sets the convergence tolerance
|
|
|
209 |
- -qep_max_it <maxits> - Sets the maximum number of iterations allowed
|
|
|
210 |
|
|
|
211 |
Notes:
|
|
|
212 |
Use PETSC_IGNORE for an argument that need not be changed.
|
|
|
213 |
|
|
|
214 |
Use PETSC_DECIDE for maxits to assign a reasonably good value, which is
|
|
|
215 |
dependent on the solution method.
|
|
|
216 |
|
|
|
217 |
Level: intermediate
|
|
|
218 |
|
|
|
219 |
.seealso: QEPGetTolerances()
|
|
|
220 |
@*/
|
|
|
221 |
PetscErrorCode QEPSetTolerances(QEP qep,PetscReal tol,PetscInt maxits)
|
|
|
222 |
{
|
|
|
223 |
PetscFunctionBegin;
|
| 2213 |
jroman |
224 |
PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
|
| 2326 |
jroman |
225 |
PetscValidLogicalCollectiveReal(qep,tol,2);
|
|
|
226 |
PetscValidLogicalCollectiveInt(qep,maxits,3);
|
| 1887 |
jroman |
227 |
if (tol != PETSC_IGNORE) {
|
|
|
228 |
if (tol == PETSC_DEFAULT) {
|
|
|
229 |
qep->tol = 1e-7;
|
|
|
230 |
} else {
|
| 2214 |
jroman |
231 |
if (tol < 0.0) SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
|
| 1887 |
jroman |
232 |
qep->tol = tol;
|
|
|
233 |
}
|
|
|
234 |
}
|
|
|
235 |
if (maxits != PETSC_IGNORE) {
|
|
|
236 |
if (maxits == PETSC_DEFAULT || maxits == PETSC_DECIDE) {
|
|
|
237 |
qep->max_it = 0;
|
|
|
238 |
qep->setupcalled = 0;
|
|
|
239 |
} else {
|
| 2214 |
jroman |
240 |
if (maxits < 0) SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of maxits. Must be > 0");
|
| 1887 |
jroman |
241 |
qep->max_it = maxits;
|
|
|
242 |
}
|
|
|
243 |
}
|
|
|
244 |
PetscFunctionReturn(0);
|
|
|
245 |
}
|
|
|
246 |
|
|
|
247 |
#undef __FUNCT__
|
|
|
248 |
#define __FUNCT__ "QEPGetDimensions"
|
|
|
249 |
/*@
|
|
|
250 |
QEPGetDimensions - Gets the number of eigenvalues to compute
|
|
|
251 |
and the dimension of the subspace.
|
|
|
252 |
|
|
|
253 |
Not Collective
|
|
|
254 |
|
|
|
255 |
Input Parameter:
|
|
|
256 |
. qep - the quadratic eigensolver context
|
|
|
257 |
|
|
|
258 |
Output Parameters:
|
|
|
259 |
+ nev - number of eigenvalues to compute
|
|
|
260 |
. ncv - the maximum dimension of the subspace to be used by the solver
|
|
|
261 |
- mpd - the maximum dimension allowed for the projected problem
|
|
|
262 |
|
|
|
263 |
Notes:
|
|
|
264 |
The user can specify PETSC_NULL for any parameter that is not needed.
|
|
|
265 |
|
|
|
266 |
Level: intermediate
|
|
|
267 |
|
|
|
268 |
.seealso: QEPSetDimensions()
|
|
|
269 |
@*/
|
|
|
270 |
PetscErrorCode QEPGetDimensions(QEP qep,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)
|
|
|
271 |
{
|
|
|
272 |
PetscFunctionBegin;
|
| 2213 |
jroman |
273 |
PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
|
| 1887 |
jroman |
274 |
if (nev) *nev = qep->nev;
|
|
|
275 |
if (ncv) *ncv = qep->ncv;
|
|
|
276 |
if (mpd) *mpd = qep->mpd;
|
|
|
277 |
PetscFunctionReturn(0);
|
|
|
278 |
}
|
|
|
279 |
|
|
|
280 |
#undef __FUNCT__
|
|
|
281 |
#define __FUNCT__ "QEPSetDimensions"
|
|
|
282 |
/*@
|
|
|
283 |
QEPSetDimensions - Sets the number of eigenvalues to compute
|
|
|
284 |
and the dimension of the subspace.
|
|
|
285 |
|
| 2328 |
jroman |
286 |
Logically Collective on QEP
|
| 1887 |
jroman |
287 |
|
|
|
288 |
Input Parameters:
|
|
|
289 |
+ qep - the quadratic eigensolver context
|
|
|
290 |
. nev - number of eigenvalues to compute
|
|
|
291 |
. ncv - the maximum dimension of the subspace to be used by the solver
|
|
|
292 |
- mpd - the maximum dimension allowed for the projected problem
|
|
|
293 |
|
|
|
294 |
Options Database Keys:
|
|
|
295 |
+ -qep_nev <nev> - Sets the number of eigenvalues
|
|
|
296 |
. -qep_ncv <ncv> - Sets the dimension of the subspace
|
|
|
297 |
- -qep_mpd <mpd> - Sets the maximum projected dimension
|
|
|
298 |
|
|
|
299 |
Notes:
|
|
|
300 |
Use PETSC_IGNORE to retain the previous value of any parameter.
|
|
|
301 |
|
|
|
302 |
Use PETSC_DECIDE for ncv and mpd to assign a reasonably good value, which is
|
|
|
303 |
dependent on the solution method.
|
|
|
304 |
|
|
|
305 |
The parameters ncv and mpd are intimately related, so that the user is advised
|
| 2242 |
jroman |
306 |
to set one of them at most. Normal usage is the following:
|
|
|
307 |
(a) In cases where nev is small, the user sets ncv (a reasonable default is 2*nev).
|
|
|
308 |
(b) In cases where nev is large, the user sets mpd.
|
| 1887 |
jroman |
309 |
|
|
|
310 |
The value of ncv should always be between nev and (nev+mpd), typically
|
|
|
311 |
ncv=nev+mpd. If nev is not too large, mpd=nev is a reasonable choice, otherwise
|
|
|
312 |
a smaller value should be used.
|
|
|
313 |
|
|
|
314 |
Level: intermediate
|
|
|
315 |
|
|
|
316 |
.seealso: QEPGetDimensions()
|
|
|
317 |
@*/
|
|
|
318 |
PetscErrorCode QEPSetDimensions(QEP qep,PetscInt nev,PetscInt ncv,PetscInt mpd)
|
|
|
319 |
{
|
|
|
320 |
PetscFunctionBegin;
|
| 2213 |
jroman |
321 |
PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
|
| 2326 |
jroman |
322 |
PetscValidLogicalCollectiveInt(qep,nev,2);
|
|
|
323 |
PetscValidLogicalCollectiveInt(qep,ncv,3);
|
|
|
324 |
PetscValidLogicalCollectiveInt(qep,mpd,4);
|
| 2331 |
jroman |
325 |
if (nev != PETSC_IGNORE) {
|
| 2214 |
jroman |
326 |
if (nev<1) SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of nev. Must be > 0");
|
| 1887 |
jroman |
327 |
qep->nev = nev;
|
|
|
328 |
qep->setupcalled = 0;
|
|
|
329 |
}
|
| 2331 |
jroman |
330 |
if (ncv != PETSC_IGNORE) {
|
| 1887 |
jroman |
331 |
if (ncv == PETSC_DECIDE || ncv == PETSC_DEFAULT) {
|
|
|
332 |
qep->ncv = 0;
|
|
|
333 |
} else {
|
| 2214 |
jroman |
334 |
if (ncv<1) SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of ncv. Must be > 0");
|
| 1887 |
jroman |
335 |
qep->ncv = ncv;
|
|
|
336 |
}
|
|
|
337 |
qep->setupcalled = 0;
|
|
|
338 |
}
|
| 2331 |
jroman |
339 |
if (mpd != PETSC_IGNORE) {
|
| 1887 |
jroman |
340 |
if (mpd == PETSC_DECIDE || mpd == PETSC_DEFAULT) {
|
|
|
341 |
qep->mpd = 0;
|
|
|
342 |
} else {
|
| 2214 |
jroman |
343 |
if (mpd<1) SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of mpd. Must be > 0");
|
| 1887 |
jroman |
344 |
qep->mpd = mpd;
|
|
|
345 |
}
|
|
|
346 |
}
|
|
|
347 |
PetscFunctionReturn(0);
|
|
|
348 |
}
|
|
|
349 |
|
|
|
350 |
#undef __FUNCT__
|
|
|
351 |
#define __FUNCT__ "QEPSetWhichEigenpairs"
|
|
|
352 |
/*@
|
|
|
353 |
QEPSetWhichEigenpairs - Specifies which portion of the spectrum is
|
|
|
354 |
to be sought.
|
|
|
355 |
|
| 2328 |
jroman |
356 |
Logically Collective on QEP
|
| 1887 |
jroman |
357 |
|
|
|
358 |
Input Parameters:
|
| 1942 |
jroman |
359 |
+ qep - eigensolver context obtained from QEPCreate()
|
| 1887 |
jroman |
360 |
- which - the portion of the spectrum to be sought
|
|
|
361 |
|
|
|
362 |
Possible values:
|
|
|
363 |
The parameter 'which' can have one of these values
|
|
|
364 |
|
|
|
365 |
+ QEP_LARGEST_MAGNITUDE - largest eigenvalues in magnitude (default)
|
|
|
366 |
. QEP_SMALLEST_MAGNITUDE - smallest eigenvalues in magnitude
|
|
|
367 |
. QEP_LARGEST_REAL - largest real parts
|
|
|
368 |
. QEP_SMALLEST_REAL - smallest real parts
|
|
|
369 |
. QEP_LARGEST_IMAGINARY - largest imaginary parts
|
|
|
370 |
- QEP_SMALLEST_IMAGINARY - smallest imaginary parts
|
|
|
371 |
|
|
|
372 |
Options Database Keys:
|
|
|
373 |
+ -qep_largest_magnitude - Sets largest eigenvalues in magnitude
|
|
|
374 |
. -qep_smallest_magnitude - Sets smallest eigenvalues in magnitude
|
|
|
375 |
. -qep_largest_real - Sets largest real parts
|
|
|
376 |
. -qep_smallest_real - Sets smallest real parts
|
|
|
377 |
. -qep_largest_imaginary - Sets largest imaginary parts
|
|
|
378 |
- -qep_smallest_imaginary - Sets smallest imaginary parts
|
|
|
379 |
|
|
|
380 |
Notes:
|
|
|
381 |
Not all eigensolvers implemented in QEP account for all the possible values
|
|
|
382 |
stated above. If SLEPc is compiled for real numbers QEP_LARGEST_IMAGINARY
|
|
|
383 |
and QEP_SMALLEST_IMAGINARY use the absolute value of the imaginary part
|
|
|
384 |
for eigenvalue selection.
|
|
|
385 |
|
|
|
386 |
Level: intermediate
|
|
|
387 |
|
|
|
388 |
.seealso: QEPGetWhichEigenpairs(), QEPWhich
|
|
|
389 |
@*/
|
|
|
390 |
PetscErrorCode QEPSetWhichEigenpairs(QEP qep,QEPWhich which)
|
|
|
391 |
{
|
|
|
392 |
PetscFunctionBegin;
|
| 2213 |
jroman |
393 |
PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
|
| 2326 |
jroman |
394 |
PetscValidLogicalCollectiveEnum(qep,which,2);
|
| 1942 |
jroman |
395 |
if (which!=PETSC_IGNORE) {
|
|
|
396 |
if (which==PETSC_DECIDE || which==PETSC_DEFAULT) qep->which = (QEPWhich)0;
|
|
|
397 |
else switch (which) {
|
|
|
398 |
case QEP_LARGEST_MAGNITUDE:
|
|
|
399 |
case QEP_SMALLEST_MAGNITUDE:
|
|
|
400 |
case QEP_LARGEST_REAL:
|
|
|
401 |
case QEP_SMALLEST_REAL:
|
|
|
402 |
case QEP_LARGEST_IMAGINARY:
|
|
|
403 |
case QEP_SMALLEST_IMAGINARY:
|
|
|
404 |
if (qep->which != which) {
|
|
|
405 |
qep->setupcalled = 0;
|
|
|
406 |
qep->which = which;
|
|
|
407 |
}
|
|
|
408 |
break;
|
|
|
409 |
default:
|
| 2214 |
jroman |
410 |
SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'which' value");
|
| 1942 |
jroman |
411 |
}
|
| 1887 |
jroman |
412 |
}
|
|
|
413 |
PetscFunctionReturn(0);
|
|
|
414 |
}
|
|
|
415 |
|
|
|
416 |
#undef __FUNCT__
|
|
|
417 |
#define __FUNCT__ "QEPGetWhichEigenpairs"
|
|
|
418 |
/*@C
|
|
|
419 |
QEPGetWhichEigenpairs - Returns which portion of the spectrum is to be
|
|
|
420 |
sought.
|
|
|
421 |
|
|
|
422 |
Not Collective
|
|
|
423 |
|
|
|
424 |
Input Parameter:
|
|
|
425 |
. qep - eigensolver context obtained from QEPCreate()
|
|
|
426 |
|
|
|
427 |
Output Parameter:
|
|
|
428 |
. which - the portion of the spectrum to be sought
|
|
|
429 |
|
|
|
430 |
Notes:
|
|
|
431 |
See QEPSetWhichEigenpairs() for possible values of 'which'.
|
|
|
432 |
|
|
|
433 |
Level: intermediate
|
|
|
434 |
|
|
|
435 |
.seealso: QEPSetWhichEigenpairs(), QEPWhich
|
|
|
436 |
@*/
|
|
|
437 |
PetscErrorCode QEPGetWhichEigenpairs(QEP qep,QEPWhich *which)
|
|
|
438 |
{
|
|
|
439 |
PetscFunctionBegin;
|
| 2213 |
jroman |
440 |
PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
|
| 1887 |
jroman |
441 |
PetscValidPointer(which,2);
|
|
|
442 |
*which = qep->which;
|
|
|
443 |
PetscFunctionReturn(0);
|
|
|
444 |
}
|
|
|
445 |
|
|
|
446 |
#undef __FUNCT__
|
| 1944 |
jroman |
447 |
#define __FUNCT__ "QEPSetLeftVectorsWanted"
|
|
|
448 |
/*@
|
|
|
449 |
QEPSetLeftVectorsWanted - Specifies which eigenvectors are required.
|
|
|
450 |
|
| 2328 |
jroman |
451 |
Logically Collective on QEP
|
| 1944 |
jroman |
452 |
|
|
|
453 |
Input Parameters:
|
|
|
454 |
+ qep - the quadratic eigensolver context
|
|
|
455 |
- leftvecs - whether left eigenvectors are required or not
|
|
|
456 |
|
|
|
457 |
Options Database Keys:
|
|
|
458 |
. -qep_left_vectors <boolean> - Sets/resets the boolean flag 'leftvecs'
|
|
|
459 |
|
|
|
460 |
Notes:
|
|
|
461 |
If the user sets leftvecs=PETSC_TRUE then the solver uses a variant of
|
|
|
462 |
the algorithm that computes both right and left eigenvectors. This is
|
|
|
463 |
usually much more costly. This option is not available in all solvers.
|
|
|
464 |
|
|
|
465 |
Level: intermediate
|
|
|
466 |
|
|
|
467 |
.seealso: QEPGetLeftVectorsWanted(), QEPGetEigenvectorLeft()
|
|
|
468 |
@*/
|
| 2216 |
jroman |
469 |
PetscErrorCode QEPSetLeftVectorsWanted(QEP qep,PetscBool leftvecs)
|
| 1944 |
jroman |
470 |
{
|
|
|
471 |
PetscFunctionBegin;
|
| 2213 |
jroman |
472 |
PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
|
| 2326 |
jroman |
473 |
PetscValidLogicalCollectiveBool(qep,leftvecs,2);
|
| 1949 |
jroman |
474 |
if (qep->leftvecs != leftvecs) {
|
|
|
475 |
qep->leftvecs = leftvecs;
|
|
|
476 |
qep->setupcalled = 0;
|
|
|
477 |
}
|
| 1944 |
jroman |
478 |
PetscFunctionReturn(0);
|
|
|
479 |
}
|
|
|
480 |
|
|
|
481 |
#undef __FUNCT__
|
|
|
482 |
#define __FUNCT__ "QEPGetLeftVectorsWanted"
|
|
|
483 |
/*@C
|
|
|
484 |
QEPGetLeftVectorsWanted - Returns the flag indicating whether left
|
|
|
485 |
eigenvectors are required or not.
|
|
|
486 |
|
|
|
487 |
Not Collective
|
|
|
488 |
|
|
|
489 |
Input Parameter:
|
|
|
490 |
. qep - the eigensolver context
|
|
|
491 |
|
|
|
492 |
Output Parameter:
|
|
|
493 |
. leftvecs - the returned flag
|
|
|
494 |
|
|
|
495 |
Level: intermediate
|
|
|
496 |
|
|
|
497 |
.seealso: QEPSetLeftVectorsWanted()
|
|
|
498 |
@*/
|
| 2216 |
jroman |
499 |
PetscErrorCode QEPGetLeftVectorsWanted(QEP qep,PetscBool *leftvecs)
|
| 1944 |
jroman |
500 |
{
|
|
|
501 |
PetscFunctionBegin;
|
| 2213 |
jroman |
502 |
PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
|
| 1944 |
jroman |
503 |
PetscValidPointer(leftvecs,2);
|
|
|
504 |
*leftvecs = qep->leftvecs;
|
|
|
505 |
PetscFunctionReturn(0);
|
|
|
506 |
}
|
|
|
507 |
|
|
|
508 |
#undef __FUNCT__
|
| 1918 |
jroman |
509 |
#define __FUNCT__ "QEPGetScaleFactor"
|
|
|
510 |
/*@
|
|
|
511 |
QEPGetScaleFactor - Gets the factor used for scaling the quadratic eigenproblem.
|
|
|
512 |
|
|
|
513 |
Not Collective
|
|
|
514 |
|
|
|
515 |
Input Parameter:
|
|
|
516 |
. qep - the quadratic eigensolver context
|
|
|
517 |
|
|
|
518 |
Output Parameters:
|
|
|
519 |
. alpha - the scaling factor
|
|
|
520 |
|
|
|
521 |
Notes:
|
|
|
522 |
If the user did not specify a scaling factor, then after QEPSolve() the
|
|
|
523 |
default value is returned.
|
|
|
524 |
|
|
|
525 |
Level: intermediate
|
|
|
526 |
|
|
|
527 |
.seealso: QEPSetScaleFactor(), QEPSolve()
|
|
|
528 |
@*/
|
|
|
529 |
PetscErrorCode QEPGetScaleFactor(QEP qep,PetscReal *alpha)
|
|
|
530 |
{
|
|
|
531 |
PetscFunctionBegin;
|
| 2213 |
jroman |
532 |
PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
|
| 2319 |
jroman |
533 |
PetscValidPointer(alpha,2);
|
|
|
534 |
*alpha = qep->sfactor;
|
| 1918 |
jroman |
535 |
PetscFunctionReturn(0);
|
|
|
536 |
}
|
|
|
537 |
|
|
|
538 |
#undef __FUNCT__
|
|
|
539 |
#define __FUNCT__ "QEPSetScaleFactor"
|
|
|
540 |
/*@
|
|
|
541 |
QEPSetScaleFactor - Sets the scaling factor to be used for scaling the
|
|
|
542 |
quadratic problem before attempting to solve.
|
|
|
543 |
|
| 2328 |
jroman |
544 |
Logically Collective on QEP
|
| 1918 |
jroman |
545 |
|
|
|
546 |
Input Parameters:
|
|
|
547 |
+ qep - the quadratic eigensolver context
|
|
|
548 |
- alpha - the scaling factor
|
|
|
549 |
|
|
|
550 |
Options Database Keys:
|
|
|
551 |
. -qep_scale <alpha> - Sets the scaling factor
|
|
|
552 |
|
|
|
553 |
Notes:
|
|
|
554 |
For the problem (l^2*M + l*C + K)*x = 0, the effect of scaling is to work
|
|
|
555 |
with matrices (alpha^2*M, alpha*C, K), then scale the computed eigenvalue.
|
|
|
556 |
|
|
|
557 |
The default is to scale with alpha = norm(K)/norm(M).
|
|
|
558 |
|
|
|
559 |
Level: intermediate
|
|
|
560 |
|
|
|
561 |
.seealso: QEPGetScaleFactor()
|
|
|
562 |
@*/
|
|
|
563 |
PetscErrorCode QEPSetScaleFactor(QEP qep,PetscReal alpha)
|
|
|
564 |
{
|
|
|
565 |
PetscFunctionBegin;
|
| 2213 |
jroman |
566 |
PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
|
| 2326 |
jroman |
567 |
PetscValidLogicalCollectiveReal(qep,alpha,2);
|
| 1918 |
jroman |
568 |
if (alpha != PETSC_IGNORE) {
|
|
|
569 |
if (alpha == PETSC_DEFAULT || alpha == PETSC_DECIDE) {
|
|
|
570 |
qep->sfactor = 0.0;
|
|
|
571 |
} else {
|
| 2214 |
jroman |
572 |
if (alpha < 0.0) SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of alpha. Must be > 0");
|
| 1918 |
jroman |
573 |
qep->sfactor = alpha;
|
|
|
574 |
}
|
|
|
575 |
}
|
|
|
576 |
PetscFunctionReturn(0);
|
|
|
577 |
}
|
|
|
578 |
|
|
|
579 |
#undef __FUNCT__
|
| 1906 |
jroman |
580 |
#define __FUNCT__ "QEPSetProblemType"
|
|
|
581 |
/*@
|
|
|
582 |
QEPSetProblemType - Specifies the type of the quadratic eigenvalue problem.
|
|
|
583 |
|
| 2328 |
jroman |
584 |
Logically Collective on QEP
|
| 1906 |
jroman |
585 |
|
|
|
586 |
Input Parameters:
|
|
|
587 |
+ qep - the quadratic eigensolver context
|
|
|
588 |
- type - a known type of quadratic eigenvalue problem
|
|
|
589 |
|
|
|
590 |
Options Database Keys:
|
|
|
591 |
+ -qep_general - general problem with no particular structure
|
|
|
592 |
. -qep_hermitian - problem whose coefficient matrices are Hermitian
|
|
|
593 |
- -qep_gyroscopic - problem with Hamiltonian structure
|
|
|
594 |
|
|
|
595 |
Notes:
|
|
|
596 |
Allowed values for the problem type are: general (QEP_GENERAL), Hermitian
|
|
|
597 |
(QEP_HERMITIAN), and gyroscopic (QEP_GYROSCOPIC).
|
|
|
598 |
|
|
|
599 |
This function is used to instruct SLEPc to exploit certain structure in
|
|
|
600 |
the quadratic eigenproblem. By default, no particular structure is assumed.
|
|
|
601 |
|
|
|
602 |
If the problem matrices are Hermitian (symmetric in the real case) or
|
|
|
603 |
Hermitian/skew-Hermitian then the solver can exploit this fact to perform
|
|
|
604 |
less operations or provide better stability.
|
|
|
605 |
|
|
|
606 |
Level: intermediate
|
|
|
607 |
|
|
|
608 |
.seealso: QEPSetOperators(), QEPSetType(), QEPGetProblemType(), QEPProblemType
|
|
|
609 |
@*/
|
|
|
610 |
PetscErrorCode QEPSetProblemType(QEP qep,QEPProblemType type)
|
|
|
611 |
{
|
|
|
612 |
PetscFunctionBegin;
|
| 2213 |
jroman |
613 |
PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
|
| 2326 |
jroman |
614 |
PetscValidLogicalCollectiveEnum(qep,type,2);
|
| 1906 |
jroman |
615 |
if (type!=QEP_GENERAL && type!=QEP_HERMITIAN && type!=QEP_GYROSCOPIC)
|
| 2214 |
jroman |
616 |
SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_WRONG,"Unknown eigenvalue problem type");
|
| 1906 |
jroman |
617 |
qep->problem_type = type;
|
|
|
618 |
PetscFunctionReturn(0);
|
|
|
619 |
}
|
|
|
620 |
|
|
|
621 |
#undef __FUNCT__
|
|
|
622 |
#define __FUNCT__ "QEPGetProblemType"
|
|
|
623 |
/*@C
|
|
|
624 |
QEPGetProblemType - Gets the problem type from the QEP object.
|
|
|
625 |
|
|
|
626 |
Not Collective
|
|
|
627 |
|
|
|
628 |
Input Parameter:
|
|
|
629 |
. qep - the quadratic eigensolver context
|
|
|
630 |
|
|
|
631 |
Output Parameter:
|
|
|
632 |
. type - name of QEP problem type
|
|
|
633 |
|
|
|
634 |
Level: intermediate
|
|
|
635 |
|
|
|
636 |
.seealso: QEPSetProblemType(), QEPProblemType
|
|
|
637 |
@*/
|
|
|
638 |
PetscErrorCode QEPGetProblemType(QEP qep,QEPProblemType *type)
|
|
|
639 |
{
|
|
|
640 |
PetscFunctionBegin;
|
| 2213 |
jroman |
641 |
PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
|
| 1906 |
jroman |
642 |
PetscValidPointer(type,2);
|
|
|
643 |
*type = qep->problem_type;
|
|
|
644 |
PetscFunctionReturn(0);
|
|
|
645 |
}
|
|
|
646 |
|
|
|
647 |
#undef __FUNCT__
|
| 1887 |
jroman |
648 |
#define __FUNCT__ "QEPSetConvergenceTest"
|
|
|
649 |
/*@C
|
| 2070 |
jroman |
650 |
QEPSetConvergenceTest - Sets a function to compute the error estimate used in
|
|
|
651 |
the convergence test.
|
| 1887 |
jroman |
652 |
|
| 2328 |
jroman |
653 |
Logically Collective on QEP
|
| 1887 |
jroman |
654 |
|
|
|
655 |
Input Parameters:
|
|
|
656 |
+ qep - eigensolver context obtained from QEPCreate()
|
|
|
657 |
. func - a pointer to the convergence test function
|
|
|
658 |
- ctx - a context pointer (the last parameter to the convergence test function)
|
|
|
659 |
|
|
|
660 |
Calling Sequence of func:
|
| 2070 |
jroman |
661 |
$ func(QEP qep,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal* errest,void *ctx)
|
| 1887 |
jroman |
662 |
|
|
|
663 |
+ qep - eigensolver context obtained from QEPCreate()
|
| 2070 |
jroman |
664 |
. eigr - real part of the eigenvalue
|
|
|
665 |
. eigi - imaginary part of the eigenvalue
|
|
|
666 |
. res - residual norm associated to the eigenpair
|
|
|
667 |
. errest - (output) computed error estimate
|
| 1887 |
jroman |
668 |
- ctx - optional context, as set by QEPSetConvergenceTest()
|
|
|
669 |
|
|
|
670 |
Note:
|
| 2070 |
jroman |
671 |
If the error estimate returned by the convergence test function is less than
|
|
|
672 |
the tolerance, then the eigenvalue is accepted as converged.
|
|
|
673 |
|
| 1887 |
jroman |
674 |
Level: advanced
|
|
|
675 |
|
|
|
676 |
.seealso: QEPSetTolerances()
|
|
|
677 |
@*/
|
| 2240 |
jroman |
678 |
extern PetscErrorCode QEPSetConvergenceTest(QEP qep,PetscErrorCode (*func)(QEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*),void* ctx)
|
| 1887 |
jroman |
679 |
{
|
|
|
680 |
PetscFunctionBegin;
|
| 2326 |
jroman |
681 |
PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
|
| 1887 |
jroman |
682 |
qep->conv_func = func;
|
|
|
683 |
qep->conv_ctx = ctx;
|
|
|
684 |
PetscFunctionReturn(0);
|
|
|
685 |
}
|
|
|
686 |
|
|
|
687 |
#undef __FUNCT__
|
| 2054 |
eromero |
688 |
#define __FUNCT__ "QEPSetTrackAll"
|
|
|
689 |
/*@
|
| 2317 |
jroman |
690 |
QEPSetTrackAll - Specifies if the solver must compute the residual of all
|
|
|
691 |
approximate eigenpairs or not.
|
| 2054 |
eromero |
692 |
|
| 2328 |
jroman |
693 |
Logically Collective on QEP
|
| 2054 |
eromero |
694 |
|
| 2317 |
jroman |
695 |
Input Parameters:
|
|
|
696 |
+ qep - the eigensolver context
|
|
|
697 |
- trackall - whether compute all residuals or not
|
| 2054 |
eromero |
698 |
|
| 2317 |
jroman |
699 |
Notes:
|
|
|
700 |
If the user sets trackall=PETSC_TRUE then the solver explicitly computes
|
|
|
701 |
the residual for each eigenpair approximation. Computing the residual is
|
|
|
702 |
usually an expensive operation and solvers commonly compute the associated
|
|
|
703 |
residual to the first unconverged eigenpair.
|
|
|
704 |
The options '-qep_monitor_all' and '-qep_monitor_draw_all' automatically
|
|
|
705 |
activates this option.
|
| 2054 |
eromero |
706 |
|
| 2317 |
jroman |
707 |
Level: intermediate
|
| 2054 |
eromero |
708 |
|
| 2242 |
jroman |
709 |
.seealso: QEPGetTrackAll()
|
| 2054 |
eromero |
710 |
@*/
|
| 2216 |
jroman |
711 |
PetscErrorCode QEPSetTrackAll(QEP qep,PetscBool trackall)
|
| 2054 |
eromero |
712 |
{
|
|
|
713 |
PetscFunctionBegin;
|
| 2213 |
jroman |
714 |
PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
|
| 2326 |
jroman |
715 |
PetscValidLogicalCollectiveBool(qep,trackall,2);
|
| 2054 |
eromero |
716 |
qep->trackall = trackall;
|
|
|
717 |
PetscFunctionReturn(0);
|
|
|
718 |
}
|
|
|
719 |
|
|
|
720 |
#undef __FUNCT__
|
|
|
721 |
#define __FUNCT__ "QEPGetTrackAll"
|
|
|
722 |
/*@
|
| 2317 |
jroman |
723 |
QEPGetTrackAll - Returns the flag indicating whether all residual norms must
|
|
|
724 |
be computed or not.
|
| 2054 |
eromero |
725 |
|
| 2317 |
jroman |
726 |
Not Collective
|
| 2054 |
eromero |
727 |
|
| 2317 |
jroman |
728 |
Input Parameter:
|
|
|
729 |
. qep - the eigensolver context
|
| 2054 |
eromero |
730 |
|
| 2317 |
jroman |
731 |
Output Parameter:
|
|
|
732 |
. trackall - the returned flag
|
| 2054 |
eromero |
733 |
|
| 2317 |
jroman |
734 |
Level: intermediate
|
| 2054 |
eromero |
735 |
|
| 2242 |
jroman |
736 |
.seealso: QEPSetTrackAll()
|
| 2054 |
eromero |
737 |
@*/
|
| 2216 |
jroman |
738 |
PetscErrorCode QEPGetTrackAll(QEP qep,PetscBool *trackall)
|
| 2054 |
eromero |
739 |
{
|
|
|
740 |
PetscFunctionBegin;
|
| 2213 |
jroman |
741 |
PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
|
| 2054 |
eromero |
742 |
PetscValidPointer(trackall,2);
|
|
|
743 |
*trackall = qep->trackall;
|
|
|
744 |
PetscFunctionReturn(0);
|
|
|
745 |
}
|
|
|
746 |
|
|
|
747 |
#undef __FUNCT__
|
| 1887 |
jroman |
748 |
#define __FUNCT__ "QEPSetOptionsPrefix"
|
|
|
749 |
/*@C
|
|
|
750 |
QEPSetOptionsPrefix - Sets the prefix used for searching for all
|
|
|
751 |
QEP options in the database.
|
|
|
752 |
|
| 2328 |
jroman |
753 |
Logically Collective on QEP
|
| 1887 |
jroman |
754 |
|
|
|
755 |
Input Parameters:
|
|
|
756 |
+ qep - the quadratic eigensolver context
|
|
|
757 |
- prefix - the prefix string to prepend to all QEP option requests
|
|
|
758 |
|
|
|
759 |
Notes:
|
|
|
760 |
A hyphen (-) must NOT be given at the beginning of the prefix name.
|
|
|
761 |
The first character of all runtime options is AUTOMATICALLY the
|
|
|
762 |
hyphen.
|
|
|
763 |
|
|
|
764 |
For example, to distinguish between the runtime options for two
|
|
|
765 |
different QEP contexts, one could call
|
|
|
766 |
.vb
|
|
|
767 |
QEPSetOptionsPrefix(qep1,"qeig1_")
|
|
|
768 |
QEPSetOptionsPrefix(qep2,"qeig2_")
|
|
|
769 |
.ve
|
|
|
770 |
|
|
|
771 |
Level: advanced
|
|
|
772 |
|
|
|
773 |
.seealso: QEPAppendOptionsPrefix(), QEPGetOptionsPrefix()
|
|
|
774 |
@*/
|
|
|
775 |
PetscErrorCode QEPSetOptionsPrefix(QEP qep,const char *prefix)
|
|
|
776 |
{
|
|
|
777 |
PetscErrorCode ierr;
|
| 2317 |
jroman |
778 |
|
| 1887 |
jroman |
779 |
PetscFunctionBegin;
|
| 2213 |
jroman |
780 |
PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
|
| 1887 |
jroman |
781 |
ierr = PetscObjectSetOptionsPrefix((PetscObject)qep,prefix);CHKERRQ(ierr);
|
|
|
782 |
ierr = IPSetOptionsPrefix(qep->ip,prefix);CHKERRQ(ierr);
|
|
|
783 |
ierr = IPAppendOptionsPrefix(qep->ip,"qep_");CHKERRQ(ierr);
|
|
|
784 |
PetscFunctionReturn(0);
|
|
|
785 |
}
|
|
|
786 |
|
|
|
787 |
#undef __FUNCT__
|
|
|
788 |
#define __FUNCT__ "QEPAppendOptionsPrefix"
|
|
|
789 |
/*@C
|
|
|
790 |
QEPAppendOptionsPrefix - Appends to the prefix used for searching for all
|
|
|
791 |
QEP options in the database.
|
|
|
792 |
|
| 2328 |
jroman |
793 |
Logically Collective on QEP
|
| 1887 |
jroman |
794 |
|
|
|
795 |
Input Parameters:
|
|
|
796 |
+ qep - the quadratic eigensolver context
|
|
|
797 |
- prefix - the prefix string to prepend to all QEP option requests
|
|
|
798 |
|
|
|
799 |
Notes:
|
|
|
800 |
A hyphen (-) must NOT be given at the beginning of the prefix name.
|
|
|
801 |
The first character of all runtime options is AUTOMATICALLY the hyphen.
|
|
|
802 |
|
|
|
803 |
Level: advanced
|
|
|
804 |
|
|
|
805 |
.seealso: QEPSetOptionsPrefix(), QEPGetOptionsPrefix()
|
|
|
806 |
@*/
|
|
|
807 |
PetscErrorCode QEPAppendOptionsPrefix(QEP qep,const char *prefix)
|
|
|
808 |
{
|
|
|
809 |
PetscErrorCode ierr;
|
| 2317 |
jroman |
810 |
|
| 1887 |
jroman |
811 |
PetscFunctionBegin;
|
| 2213 |
jroman |
812 |
PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
|
| 2331 |
jroman |
813 |
ierr = PetscObjectAppendOptionsPrefix((PetscObject)qep,prefix);CHKERRQ(ierr);
|
| 1887 |
jroman |
814 |
ierr = IPSetOptionsPrefix(qep->ip,prefix);CHKERRQ(ierr);
|
|
|
815 |
ierr = IPAppendOptionsPrefix(qep->ip,"qep_");CHKERRQ(ierr);
|
|
|
816 |
PetscFunctionReturn(0);
|
|
|
817 |
}
|
|
|
818 |
|
|
|
819 |
#undef __FUNCT__
|
|
|
820 |
#define __FUNCT__ "QEPGetOptionsPrefix"
|
|
|
821 |
/*@C
|
|
|
822 |
QEPGetOptionsPrefix - Gets the prefix used for searching for all
|
|
|
823 |
QEP options in the database.
|
|
|
824 |
|
|
|
825 |
Not Collective
|
|
|
826 |
|
|
|
827 |
Input Parameters:
|
|
|
828 |
. qep - the quadratic eigensolver context
|
|
|
829 |
|
|
|
830 |
Output Parameters:
|
|
|
831 |
. prefix - pointer to the prefix string used is returned
|
|
|
832 |
|
|
|
833 |
Notes: On the fortran side, the user should pass in a string 'prefix' of
|
|
|
834 |
sufficient length to hold the prefix.
|
|
|
835 |
|
|
|
836 |
Level: advanced
|
|
|
837 |
|
|
|
838 |
.seealso: QEPSetOptionsPrefix(), QEPAppendOptionsPrefix()
|
|
|
839 |
@*/
|
|
|
840 |
PetscErrorCode QEPGetOptionsPrefix(QEP qep,const char *prefix[])
|
|
|
841 |
{
|
|
|
842 |
PetscErrorCode ierr;
|
| 2317 |
jroman |
843 |
|
| 1887 |
jroman |
844 |
PetscFunctionBegin;
|
| 2213 |
jroman |
845 |
PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
|
| 1887 |
jroman |
846 |
PetscValidPointer(prefix,2);
|
|
|
847 |
ierr = PetscObjectGetOptionsPrefix((PetscObject)qep,prefix);CHKERRQ(ierr);
|
|
|
848 |
PetscFunctionReturn(0);
|
|
|
849 |
}
|