Line data Source code
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 : NEP routines related to various viewers
12 : */
13 :
14 : #include <slepc/private/nepimpl.h> /*I "slepcnep.h" I*/
15 : #include <slepc/private/bvimpl.h>
16 : #include <petscdraw.h>
17 :
18 : /*@
19 : NEPView - Prints the NEP data structure.
20 :
21 : Collective
22 :
23 : Input Parameters:
24 : + nep - the nonlinear eigenproblem solver context
25 : - viewer - optional visualization context
26 :
27 : Options Database Key:
28 : . -nep_view - Calls NEPView() at end of NEPSolve()
29 :
30 : Note:
31 : The available visualization contexts include
32 : + PETSC_VIEWER_STDOUT_SELF - standard output (default)
33 : - PETSC_VIEWER_STDOUT_WORLD - synchronized standard
34 : output where only the first processor opens
35 : the file. All other processors send their
36 : data to the first processor to print.
37 :
38 : The user can open an alternative visualization context with
39 : PetscViewerASCIIOpen() - output to a specified file.
40 :
41 : Level: beginner
42 :
43 : .seealso: FNView()
44 : @*/
45 5 : PetscErrorCode NEPView(NEP nep,PetscViewer viewer)
46 : {
47 5 : const char *type=NULL;
48 5 : char str[50];
49 5 : PetscInt i;
50 5 : PetscBool isascii,istrivial;
51 5 : PetscViewer sviewer;
52 5 : MPI_Comm child;
53 :
54 5 : PetscFunctionBegin;
55 5 : PetscValidHeaderSpecific(nep,NEP_CLASSID,1);
56 5 : if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)nep),&viewer));
57 5 : PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
58 5 : PetscCheckSameComm(nep,1,viewer,2);
59 :
60 5 : PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
61 5 : if (isascii) {
62 5 : PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)nep,viewer));
63 5 : PetscCall(PetscViewerASCIIPushTab(viewer));
64 5 : PetscTryTypeMethod(nep,view,viewer);
65 5 : PetscCall(PetscViewerASCIIPopTab(viewer));
66 5 : if (nep->problem_type) {
67 5 : switch (nep->problem_type) {
68 4 : case NEP_GENERAL: type = "general nonlinear eigenvalue problem"; break;
69 1 : case NEP_RATIONAL: type = "rational eigenvalue problem"; break;
70 : }
71 : } else type = "not yet set";
72 5 : PetscCall(PetscViewerASCIIPrintf(viewer," problem type: %s\n",type));
73 5 : if (nep->fui) {
74 5 : switch (nep->fui) {
75 0 : case NEP_USER_INTERFACE_CALLBACK:
76 0 : PetscCall(PetscViewerASCIIPrintf(viewer," nonlinear operator from user callbacks\n"));
77 : break;
78 5 : case NEP_USER_INTERFACE_SPLIT:
79 5 : PetscCall(PetscViewerASCIIPrintf(viewer," nonlinear operator in split form\n"));
80 5 : PetscCall(PetscViewerASCIIPrintf(viewer," number of terms: %" PetscInt_FMT "\n",nep->nt));
81 5 : PetscCall(PetscViewerASCIIPrintf(viewer," nonzero pattern of the matrices: %s\n",MatStructures[nep->mstr]));
82 : break;
83 : }
84 0 : } else PetscCall(PetscViewerASCIIPrintf(viewer," nonlinear operator not specified yet\n"));
85 5 : PetscCall(PetscViewerASCIIPrintf(viewer," selected portion of the spectrum: "));
86 5 : PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
87 5 : PetscCall(SlepcSNPrintfScalar(str,sizeof(str),nep->target,PETSC_FALSE));
88 5 : if (!nep->which) PetscCall(PetscViewerASCIIPrintf(viewer,"not yet set\n"));
89 5 : else switch (nep->which) {
90 0 : case NEP_WHICH_USER:
91 0 : PetscCall(PetscViewerASCIIPrintf(viewer,"user defined\n"));
92 : break;
93 4 : case NEP_TARGET_MAGNITUDE:
94 4 : PetscCall(PetscViewerASCIIPrintf(viewer,"closest to target: %s (in magnitude)\n",str));
95 : break;
96 0 : case NEP_TARGET_REAL:
97 0 : PetscCall(PetscViewerASCIIPrintf(viewer,"closest to target: %s (along the real axis)\n",str));
98 : break;
99 0 : case NEP_TARGET_IMAGINARY:
100 0 : PetscCall(PetscViewerASCIIPrintf(viewer,"closest to target: %s (along the imaginary axis)\n",str));
101 : break;
102 0 : case NEP_LARGEST_MAGNITUDE:
103 0 : PetscCall(PetscViewerASCIIPrintf(viewer,"largest eigenvalues in magnitude\n"));
104 : break;
105 0 : case NEP_SMALLEST_MAGNITUDE:
106 0 : PetscCall(PetscViewerASCIIPrintf(viewer,"smallest eigenvalues in magnitude\n"));
107 : break;
108 0 : case NEP_LARGEST_REAL:
109 0 : PetscCall(PetscViewerASCIIPrintf(viewer,"largest real parts\n"));
110 : break;
111 0 : case NEP_SMALLEST_REAL:
112 0 : PetscCall(PetscViewerASCIIPrintf(viewer,"smallest real parts\n"));
113 : break;
114 0 : case NEP_LARGEST_IMAGINARY:
115 0 : PetscCall(PetscViewerASCIIPrintf(viewer,"largest imaginary parts\n"));
116 : break;
117 0 : case NEP_SMALLEST_IMAGINARY:
118 0 : PetscCall(PetscViewerASCIIPrintf(viewer,"smallest imaginary parts\n"));
119 : break;
120 1 : case NEP_ALL:
121 1 : PetscCall(PetscViewerASCIIPrintf(viewer,"all eigenvalues in the region\n"));
122 : break;
123 : }
124 5 : PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
125 5 : if (nep->twosided) PetscCall(PetscViewerASCIIPrintf(viewer," using two-sided variant (for left eigenvectors)\n"));
126 5 : PetscCall(PetscViewerASCIIPrintf(viewer," number of eigenvalues (nev): %" PetscInt_FMT "\n",nep->nev));
127 5 : PetscCall(PetscViewerASCIIPrintf(viewer," number of column vectors (ncv): %" PetscInt_FMT "\n",nep->ncv));
128 5 : PetscCall(PetscViewerASCIIPrintf(viewer," maximum dimension of projected problem (mpd): %" PetscInt_FMT "\n",nep->mpd));
129 5 : PetscCall(PetscViewerASCIIPrintf(viewer," maximum number of iterations: %" PetscInt_FMT "\n",nep->max_it));
130 5 : PetscCall(PetscViewerASCIIPrintf(viewer," tolerance: %g\n",(double)nep->tol));
131 5 : PetscCall(PetscViewerASCIIPrintf(viewer," convergence test: "));
132 5 : PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
133 5 : switch (nep->conv) {
134 1 : case NEP_CONV_ABS:
135 1 : PetscCall(PetscViewerASCIIPrintf(viewer,"absolute\n"));break;
136 4 : case NEP_CONV_REL:
137 4 : PetscCall(PetscViewerASCIIPrintf(viewer,"relative to the eigenvalue\n"));break;
138 0 : case NEP_CONV_NORM:
139 0 : PetscCall(PetscViewerASCIIPrintf(viewer,"relative to the matrix norms\n"));
140 0 : if (nep->nrma) {
141 0 : PetscCall(PetscViewerASCIIPrintf(viewer," computed matrix norms: %g",(double)nep->nrma[0]));
142 0 : for (i=1;i<nep->nt;i++) PetscCall(PetscViewerASCIIPrintf(viewer,", %g",(double)nep->nrma[i]));
143 0 : PetscCall(PetscViewerASCIIPrintf(viewer,"\n"));
144 : }
145 : break;
146 0 : case NEP_CONV_USER:
147 0 : PetscCall(PetscViewerASCIIPrintf(viewer,"user-defined\n"));break;
148 : }
149 5 : PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
150 5 : if (nep->refine) {
151 1 : PetscCall(PetscViewerASCIIPrintf(viewer," iterative refinement: %s, with %s scheme\n",NEPRefineTypes[nep->refine],NEPRefineSchemes[nep->scheme]));
152 1 : PetscCall(PetscViewerASCIIPrintf(viewer," refinement stopping criterion: tol=%g, its=%" PetscInt_FMT "\n",(double)nep->rtol,nep->rits));
153 1 : if (nep->npart>1) PetscCall(PetscViewerASCIIPrintf(viewer," splitting communicator in %" PetscInt_FMT " partitions for refinement\n",nep->npart));
154 : }
155 5 : if (nep->nini) PetscCall(PetscViewerASCIIPrintf(viewer," dimension of user-provided initial space: %" PetscInt_FMT "\n",PetscAbs(nep->nini)));
156 0 : } else PetscTryTypeMethod(nep,view,viewer);
157 5 : PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO));
158 5 : if (!nep->V) PetscCall(NEPGetBV(nep,&nep->V));
159 5 : PetscCall(BVView(nep->V,viewer));
160 5 : if (!nep->rg) PetscCall(NEPGetRG(nep,&nep->rg));
161 5 : PetscCall(RGIsTrivial(nep->rg,&istrivial));
162 5 : if (!istrivial) PetscCall(RGView(nep->rg,viewer));
163 5 : if (nep->useds) {
164 5 : if (!nep->ds) PetscCall(NEPGetDS(nep,&nep->ds));
165 5 : PetscCall(DSView(nep->ds,viewer));
166 : }
167 5 : PetscCall(PetscViewerPopFormat(viewer));
168 5 : if (nep->refine!=NEP_REFINE_NONE) {
169 1 : PetscCall(PetscViewerASCIIPushTab(viewer));
170 1 : if (nep->npart>1) {
171 0 : if (nep->refinesubc->color==0) {
172 0 : PetscCall(PetscSubcommGetChild(nep->refinesubc,&child));
173 0 : PetscCall(PetscViewerASCIIGetStdout(child,&sviewer));
174 0 : PetscCall(KSPView(nep->refineksp,sviewer));
175 : }
176 1 : } else PetscCall(KSPView(nep->refineksp,viewer));
177 1 : PetscCall(PetscViewerASCIIPopTab(viewer));
178 : }
179 5 : PetscFunctionReturn(PETSC_SUCCESS);
180 : }
181 :
182 : /*@
183 : NEPViewFromOptions - View from options
184 :
185 : Collective
186 :
187 : Input Parameters:
188 : + nep - the nonlinear eigensolver context
189 : . obj - optional object
190 : - name - command line option
191 :
192 : Level: intermediate
193 :
194 : .seealso: NEPView(), NEPCreate()
195 : @*/
196 320 : PetscErrorCode NEPViewFromOptions(NEP nep,PetscObject obj,const char name[])
197 : {
198 320 : PetscFunctionBegin;
199 320 : PetscValidHeaderSpecific(nep,NEP_CLASSID,1);
200 320 : PetscCall(PetscObjectViewFromOptions((PetscObject)nep,obj,name));
201 320 : PetscFunctionReturn(PETSC_SUCCESS);
202 : }
203 :
204 : /*@
205 : NEPConvergedReasonView - Displays the reason a NEP solve converged or diverged.
206 :
207 : Collective
208 :
209 : Input Parameters:
210 : + nep - the nonlinear eigensolver context
211 : - viewer - the viewer to display the reason
212 :
213 : Options Database Keys:
214 : . -nep_converged_reason - print reason for convergence, and number of iterations
215 :
216 : Note:
217 : To change the format of the output call PetscViewerPushFormat(viewer,format) before
218 : this call. Use PETSC_VIEWER_DEFAULT for the default, use PETSC_VIEWER_FAILED to only
219 : display a reason if it fails. The latter can be set in the command line with
220 : -nep_converged_reason ::failed
221 :
222 : Level: intermediate
223 :
224 : .seealso: NEPSetConvergenceTest(), NEPSetTolerances(), NEPGetIterationNumber(), NEPConvergedReasonViewFromOptions()
225 : @*/
226 2 : PetscErrorCode NEPConvergedReasonView(NEP nep,PetscViewer viewer)
227 : {
228 2 : PetscBool isAscii;
229 2 : PetscViewerFormat format;
230 :
231 2 : PetscFunctionBegin;
232 2 : if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)nep));
233 2 : PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isAscii));
234 2 : if (isAscii) {
235 2 : PetscCall(PetscViewerGetFormat(viewer,&format));
236 2 : PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)nep)->tablevel));
237 3 : if (nep->reason > 0 && format != PETSC_VIEWER_FAILED) PetscCall(PetscViewerASCIIPrintf(viewer,"%s Nonlinear eigensolve converged (%" PetscInt_FMT " eigenpair%s) due to %s; iterations %" PetscInt_FMT "\n",((PetscObject)nep)->prefix?((PetscObject)nep)->prefix:"",nep->nconv,(nep->nconv>1)?"s":"",NEPConvergedReasons[nep->reason],nep->its));
238 0 : else if (nep->reason <= 0) PetscCall(PetscViewerASCIIPrintf(viewer,"%s Nonlinear eigensolve did not converge due to %s; iterations %" PetscInt_FMT "\n",((PetscObject)nep)->prefix?((PetscObject)nep)->prefix:"",NEPConvergedReasons[nep->reason],nep->its));
239 2 : PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)nep)->tablevel));
240 : }
241 2 : PetscFunctionReturn(PETSC_SUCCESS);
242 : }
243 :
244 : /*@
245 : NEPConvergedReasonViewFromOptions - Processes command line options to determine if/how
246 : the NEP converged reason is to be viewed.
247 :
248 : Collective
249 :
250 : Input Parameter:
251 : . nep - the nonlinear eigensolver context
252 :
253 : Level: developer
254 :
255 : .seealso: NEPConvergedReasonView()
256 : @*/
257 160 : PetscErrorCode NEPConvergedReasonViewFromOptions(NEP nep)
258 : {
259 160 : PetscViewer viewer;
260 160 : PetscBool flg;
261 160 : static PetscBool incall = PETSC_FALSE;
262 160 : PetscViewerFormat format;
263 :
264 160 : PetscFunctionBegin;
265 160 : if (incall) PetscFunctionReturn(PETSC_SUCCESS);
266 160 : incall = PETSC_TRUE;
267 160 : PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)nep),((PetscObject)nep)->options,((PetscObject)nep)->prefix,"-nep_converged_reason",&viewer,&format,&flg));
268 160 : if (flg) {
269 1 : PetscCall(PetscViewerPushFormat(viewer,format));
270 1 : PetscCall(NEPConvergedReasonView(nep,viewer));
271 1 : PetscCall(PetscViewerPopFormat(viewer));
272 1 : PetscCall(PetscViewerDestroy(&viewer));
273 : }
274 160 : incall = PETSC_FALSE;
275 160 : PetscFunctionReturn(PETSC_SUCCESS);
276 : }
277 :
278 153 : static PetscErrorCode NEPErrorView_ASCII(NEP nep,NEPErrorType etype,PetscViewer viewer)
279 : {
280 153 : PetscReal error;
281 153 : PetscInt i,j,k,nvals;
282 :
283 153 : PetscFunctionBegin;
284 153 : nvals = (nep->which==NEP_ALL)? nep->nconv: nep->nev;
285 153 : if (nep->which!=NEP_ALL && nep->nconv<nvals) {
286 0 : PetscCall(PetscViewerASCIIPrintf(viewer," Problem: less than %" PetscInt_FMT " eigenvalues converged\n\n",nep->nev));
287 0 : PetscFunctionReturn(PETSC_SUCCESS);
288 : }
289 153 : if (nep->which==NEP_ALL && !nvals) {
290 0 : PetscCall(PetscViewerASCIIPrintf(viewer," No eigenvalues have been found\n\n"));
291 0 : PetscFunctionReturn(PETSC_SUCCESS);
292 : }
293 622 : for (i=0;i<nvals;i++) {
294 469 : PetscCall(NEPComputeError(nep,i,etype,&error));
295 469 : if (error>=5.0*nep->tol) {
296 0 : PetscCall(PetscViewerASCIIPrintf(viewer," Problem: some of the first %" PetscInt_FMT " relative errors are higher than the tolerance\n\n",nvals));
297 469 : PetscFunctionReturn(PETSC_SUCCESS);
298 : }
299 : }
300 153 : if (nep->which==NEP_ALL) PetscCall(PetscViewerASCIIPrintf(viewer," Found %" PetscInt_FMT " eigenvalues, all of them computed up to the required tolerance:",nvals));
301 113 : else PetscCall(PetscViewerASCIIPrintf(viewer," All requested eigenvalues computed up to the required tolerance:"));
302 309 : for (i=0;i<=(nvals-1)/8;i++) {
303 156 : PetscCall(PetscViewerASCIIPrintf(viewer,"\n "));
304 625 : for (j=0;j<PetscMin(8,nvals-8*i);j++) {
305 469 : k = nep->perm[8*i+j];
306 469 : PetscCall(SlepcPrintEigenvalueASCII(viewer,nep->eigr[k],nep->eigi[k]));
307 469 : if (8*i+j+1<nvals) PetscCall(PetscViewerASCIIPrintf(viewer,", "));
308 : }
309 : }
310 153 : PetscCall(PetscViewerASCIIPrintf(viewer,"\n\n"));
311 153 : PetscFunctionReturn(PETSC_SUCCESS);
312 : }
313 :
314 1 : static PetscErrorCode NEPErrorView_DETAIL(NEP nep,NEPErrorType etype,PetscViewer viewer)
315 : {
316 1 : PetscReal error,re,im;
317 1 : PetscScalar kr,ki;
318 1 : PetscInt i;
319 1 : char ex[30],sep[]=" ---------------------- --------------------\n";
320 :
321 1 : PetscFunctionBegin;
322 1 : if (!nep->nconv) PetscFunctionReturn(PETSC_SUCCESS);
323 1 : switch (etype) {
324 0 : case NEP_ERROR_ABSOLUTE:
325 0 : PetscCall(PetscSNPrintf(ex,sizeof(ex)," ||T(k)x||"));
326 : break;
327 0 : case NEP_ERROR_RELATIVE:
328 0 : PetscCall(PetscSNPrintf(ex,sizeof(ex)," ||T(k)x||/||kx||"));
329 : break;
330 1 : case NEP_ERROR_BACKWARD:
331 1 : PetscCall(PetscSNPrintf(ex,sizeof(ex)," eta(x,k)"));
332 : break;
333 : }
334 1 : PetscCall(PetscViewerASCIIPrintf(viewer,"%s k %s\n%s",sep,ex,sep));
335 2 : for (i=0;i<nep->nconv;i++) {
336 1 : PetscCall(NEPGetEigenpair(nep,i,&kr,&ki,NULL,NULL));
337 1 : PetscCall(NEPComputeError(nep,i,etype,&error));
338 : #if defined(PETSC_USE_COMPLEX)
339 1 : re = PetscRealPart(kr);
340 1 : im = PetscImaginaryPart(kr);
341 : #else
342 : re = kr;
343 : im = ki;
344 : #endif
345 1 : if (im!=0.0) PetscCall(PetscViewerASCIIPrintf(viewer," % 9f%+9fi %12g\n",(double)re,(double)im,(double)error));
346 1 : else PetscCall(PetscViewerASCIIPrintf(viewer," % 12f %12g\n",(double)re,(double)error));
347 : }
348 1 : PetscCall(PetscViewerASCIIPrintf(viewer,"%s",sep));
349 1 : PetscFunctionReturn(PETSC_SUCCESS);
350 : }
351 :
352 2 : static PetscErrorCode NEPErrorView_MATLAB(NEP nep,NEPErrorType etype,PetscViewer viewer)
353 : {
354 2 : PetscReal error;
355 2 : PetscInt i;
356 2 : const char *name;
357 :
358 2 : PetscFunctionBegin;
359 2 : PetscCall(PetscObjectGetName((PetscObject)nep,&name));
360 2 : PetscCall(PetscViewerASCIIPrintf(viewer,"Error_%s = [\n",name));
361 7 : for (i=0;i<nep->nconv;i++) {
362 5 : PetscCall(NEPComputeError(nep,i,etype,&error));
363 5 : PetscCall(PetscViewerASCIIPrintf(viewer,"%18.16e\n",(double)error));
364 : }
365 2 : PetscCall(PetscViewerASCIIPrintf(viewer,"];\n"));
366 2 : PetscFunctionReturn(PETSC_SUCCESS);
367 : }
368 :
369 : /*@
370 : NEPErrorView - Displays the errors associated with the computed solution
371 : (as well as the eigenvalues).
372 :
373 : Collective
374 :
375 : Input Parameters:
376 : + nep - the nonlinear eigensolver context
377 : . etype - error type
378 : - viewer - optional visualization context
379 :
380 : Options Database Keys:
381 : + -nep_error_absolute - print absolute errors of each eigenpair
382 : . -nep_error_relative - print relative errors of each eigenpair
383 : - -nep_error_backward - print backward errors of each eigenpair
384 :
385 : Notes:
386 : By default, this function checks the error of all eigenpairs and prints
387 : the eigenvalues if all of them are below the requested tolerance.
388 : If the viewer has format=PETSC_VIEWER_ASCII_INFO_DETAIL then a table with
389 : eigenvalues and corresponding errors is printed.
390 :
391 : Level: intermediate
392 :
393 : .seealso: NEPSolve(), NEPValuesView(), NEPVectorsView()
394 : @*/
395 156 : PetscErrorCode NEPErrorView(NEP nep,NEPErrorType etype,PetscViewer viewer)
396 : {
397 156 : PetscBool isascii;
398 156 : PetscViewerFormat format;
399 :
400 156 : PetscFunctionBegin;
401 156 : PetscValidHeaderSpecific(nep,NEP_CLASSID,1);
402 156 : if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)nep),&viewer));
403 156 : PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,3);
404 156 : PetscCheckSameComm(nep,1,viewer,3);
405 156 : NEPCheckSolved(nep,1);
406 156 : PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
407 156 : if (!isascii) PetscFunctionReturn(PETSC_SUCCESS);
408 :
409 156 : PetscCall(PetscViewerGetFormat(viewer,&format));
410 156 : switch (format) {
411 153 : case PETSC_VIEWER_DEFAULT:
412 : case PETSC_VIEWER_ASCII_INFO:
413 153 : PetscCall(NEPErrorView_ASCII(nep,etype,viewer));
414 : break;
415 1 : case PETSC_VIEWER_ASCII_INFO_DETAIL:
416 1 : PetscCall(NEPErrorView_DETAIL(nep,etype,viewer));
417 : break;
418 2 : case PETSC_VIEWER_ASCII_MATLAB:
419 2 : PetscCall(NEPErrorView_MATLAB(nep,etype,viewer));
420 : break;
421 0 : default:
422 0 : PetscCall(PetscInfo(nep,"Unsupported viewer format %s\n",PetscViewerFormats[format]));
423 : }
424 156 : PetscFunctionReturn(PETSC_SUCCESS);
425 : }
426 :
427 : /*@
428 : NEPErrorViewFromOptions - Processes command line options to determine if/how
429 : the errors of the computed solution are to be viewed.
430 :
431 : Collective
432 :
433 : Input Parameter:
434 : . nep - the nonlinear eigensolver context
435 :
436 : Level: developer
437 :
438 : .seealso: NEPErrorView()
439 : @*/
440 160 : PetscErrorCode NEPErrorViewFromOptions(NEP nep)
441 : {
442 160 : PetscViewer viewer;
443 160 : PetscBool flg;
444 160 : static PetscBool incall = PETSC_FALSE;
445 160 : PetscViewerFormat format;
446 :
447 160 : PetscFunctionBegin;
448 160 : if (incall) PetscFunctionReturn(PETSC_SUCCESS);
449 160 : incall = PETSC_TRUE;
450 160 : PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)nep),((PetscObject)nep)->options,((PetscObject)nep)->prefix,"-nep_error_absolute",&viewer,&format,&flg));
451 160 : if (flg) {
452 1 : PetscCall(PetscViewerPushFormat(viewer,format));
453 1 : PetscCall(NEPErrorView(nep,NEP_ERROR_ABSOLUTE,viewer));
454 1 : PetscCall(PetscViewerPopFormat(viewer));
455 1 : PetscCall(PetscViewerDestroy(&viewer));
456 : }
457 160 : PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)nep),((PetscObject)nep)->options,((PetscObject)nep)->prefix,"-nep_error_relative",&viewer,&format,&flg));
458 160 : if (flg) {
459 1 : PetscCall(PetscViewerPushFormat(viewer,format));
460 1 : PetscCall(NEPErrorView(nep,NEP_ERROR_RELATIVE,viewer));
461 1 : PetscCall(PetscViewerPopFormat(viewer));
462 1 : PetscCall(PetscViewerDestroy(&viewer));
463 : }
464 160 : PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)nep),((PetscObject)nep)->options,((PetscObject)nep)->prefix,"-nep_error_backward",&viewer,&format,&flg));
465 160 : if (flg) {
466 1 : PetscCall(PetscViewerPushFormat(viewer,format));
467 1 : PetscCall(NEPErrorView(nep,NEP_ERROR_BACKWARD,viewer));
468 1 : PetscCall(PetscViewerPopFormat(viewer));
469 1 : PetscCall(PetscViewerDestroy(&viewer));
470 : }
471 160 : incall = PETSC_FALSE;
472 160 : PetscFunctionReturn(PETSC_SUCCESS);
473 : }
474 :
475 0 : static PetscErrorCode NEPValuesView_DRAW(NEP nep,PetscViewer viewer)
476 : {
477 0 : PetscDraw draw;
478 0 : PetscDrawSP drawsp;
479 0 : PetscReal re,im;
480 0 : PetscInt i,k;
481 :
482 0 : PetscFunctionBegin;
483 0 : if (!nep->nconv) PetscFunctionReturn(PETSC_SUCCESS);
484 0 : PetscCall(PetscViewerDrawGetDraw(viewer,0,&draw));
485 0 : PetscCall(PetscDrawSetTitle(draw,"Computed Eigenvalues"));
486 0 : PetscCall(PetscDrawSPCreate(draw,1,&drawsp));
487 0 : for (i=0;i<nep->nconv;i++) {
488 0 : k = nep->perm[i];
489 : #if defined(PETSC_USE_COMPLEX)
490 0 : re = PetscRealPart(nep->eigr[k]);
491 0 : im = PetscImaginaryPart(nep->eigr[k]);
492 : #else
493 : re = nep->eigr[k];
494 : im = nep->eigi[k];
495 : #endif
496 0 : PetscCall(PetscDrawSPAddPoint(drawsp,&re,&im));
497 : }
498 0 : PetscCall(PetscDrawSPDraw(drawsp,PETSC_TRUE));
499 0 : PetscCall(PetscDrawSPSave(drawsp));
500 0 : PetscCall(PetscDrawSPDestroy(&drawsp));
501 0 : PetscFunctionReturn(PETSC_SUCCESS);
502 : }
503 :
504 1 : static PetscErrorCode NEPValuesView_BINARY(NEP nep,PetscViewer viewer)
505 : {
506 : #if defined(PETSC_HAVE_COMPLEX)
507 1 : PetscInt i,k;
508 1 : PetscComplex *ev;
509 : #endif
510 :
511 1 : PetscFunctionBegin;
512 : #if defined(PETSC_HAVE_COMPLEX)
513 1 : PetscCall(PetscMalloc1(nep->nconv,&ev));
514 5 : for (i=0;i<nep->nconv;i++) {
515 4 : k = nep->perm[i];
516 : #if defined(PETSC_USE_COMPLEX)
517 4 : ev[i] = nep->eigr[k];
518 : #else
519 : ev[i] = PetscCMPLX(nep->eigr[k],nep->eigi[k]);
520 : #endif
521 : }
522 1 : PetscCall(PetscViewerBinaryWrite(viewer,ev,nep->nconv,PETSC_COMPLEX));
523 1 : PetscCall(PetscFree(ev));
524 : #endif
525 1 : PetscFunctionReturn(PETSC_SUCCESS);
526 : }
527 :
528 : #if defined(PETSC_HAVE_HDF5)
529 : static PetscErrorCode NEPValuesView_HDF5(NEP nep,PetscViewer viewer)
530 : {
531 : PetscInt i,k,n,N;
532 : PetscMPIInt rank;
533 : Vec v;
534 : char vname[30];
535 : const char *ename;
536 :
537 : PetscFunctionBegin;
538 : PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)nep),&rank));
539 : N = nep->nconv;
540 : n = rank? 0: N;
541 : /* create a vector containing the eigenvalues */
542 : PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)nep),n,N,&v));
543 : PetscCall(PetscObjectGetName((PetscObject)nep,&ename));
544 : PetscCall(PetscSNPrintf(vname,sizeof(vname),"eigr_%s",ename));
545 : PetscCall(PetscObjectSetName((PetscObject)v,vname));
546 : if (!rank) {
547 : for (i=0;i<nep->nconv;i++) {
548 : k = nep->perm[i];
549 : PetscCall(VecSetValue(v,i,nep->eigr[k],INSERT_VALUES));
550 : }
551 : }
552 : PetscCall(VecAssemblyBegin(v));
553 : PetscCall(VecAssemblyEnd(v));
554 : PetscCall(VecView(v,viewer));
555 : #if !defined(PETSC_USE_COMPLEX)
556 : /* in real scalars write the imaginary part as a separate vector */
557 : PetscCall(PetscSNPrintf(vname,sizeof(vname),"eigi_%s",ename));
558 : PetscCall(PetscObjectSetName((PetscObject)v,vname));
559 : if (!rank) {
560 : for (i=0;i<nep->nconv;i++) {
561 : k = nep->perm[i];
562 : PetscCall(VecSetValue(v,i,nep->eigi[k],INSERT_VALUES));
563 : }
564 : }
565 : PetscCall(VecAssemblyBegin(v));
566 : PetscCall(VecAssemblyEnd(v));
567 : PetscCall(VecView(v,viewer));
568 : #endif
569 : PetscCall(VecDestroy(&v));
570 : PetscFunctionReturn(PETSC_SUCCESS);
571 : }
572 : #endif
573 :
574 1 : static PetscErrorCode NEPValuesView_ASCII(NEP nep,PetscViewer viewer)
575 : {
576 1 : PetscInt i,k;
577 :
578 1 : PetscFunctionBegin;
579 1 : PetscCall(PetscViewerASCIIPrintf(viewer,"Eigenvalues = \n"));
580 2 : for (i=0;i<nep->nconv;i++) {
581 1 : k = nep->perm[i];
582 1 : PetscCall(PetscViewerASCIIPrintf(viewer," "));
583 1 : PetscCall(SlepcPrintEigenvalueASCII(viewer,nep->eigr[k],nep->eigi[k]));
584 1 : PetscCall(PetscViewerASCIIPrintf(viewer,"\n"));
585 : }
586 1 : PetscCall(PetscViewerASCIIPrintf(viewer,"\n"));
587 1 : PetscFunctionReturn(PETSC_SUCCESS);
588 : }
589 :
590 1 : static PetscErrorCode NEPValuesView_MATLAB(NEP nep,PetscViewer viewer)
591 : {
592 1 : PetscInt i,k;
593 1 : PetscReal re,im;
594 1 : const char *name;
595 :
596 1 : PetscFunctionBegin;
597 1 : PetscCall(PetscObjectGetName((PetscObject)nep,&name));
598 1 : PetscCall(PetscViewerASCIIPrintf(viewer,"Lambda_%s = [\n",name));
599 2 : for (i=0;i<nep->nconv;i++) {
600 1 : k = nep->perm[i];
601 : #if defined(PETSC_USE_COMPLEX)
602 1 : re = PetscRealPart(nep->eigr[k]);
603 1 : im = PetscImaginaryPart(nep->eigr[k]);
604 : #else
605 : re = nep->eigr[k];
606 : im = nep->eigi[k];
607 : #endif
608 1 : if (im!=0.0) PetscCall(PetscViewerASCIIPrintf(viewer,"%18.16e%+18.16ei\n",(double)re,(double)im));
609 1 : else PetscCall(PetscViewerASCIIPrintf(viewer,"%18.16e\n",(double)re));
610 : }
611 1 : PetscCall(PetscViewerASCIIPrintf(viewer,"];\n"));
612 1 : PetscFunctionReturn(PETSC_SUCCESS);
613 : }
614 :
615 : /*@
616 : NEPValuesView - Displays the computed eigenvalues in a viewer.
617 :
618 : Collective
619 :
620 : Input Parameters:
621 : + nep - the nonlinear eigensolver context
622 : - viewer - the viewer
623 :
624 : Options Database Key:
625 : . -nep_view_values - print computed eigenvalues
626 :
627 : Level: intermediate
628 :
629 : .seealso: NEPSolve(), NEPVectorsView(), NEPErrorView()
630 : @*/
631 3 : PetscErrorCode NEPValuesView(NEP nep,PetscViewer viewer)
632 : {
633 3 : PetscBool isascii,isdraw,isbinary;
634 3 : PetscViewerFormat format;
635 : #if defined(PETSC_HAVE_HDF5)
636 : PetscBool ishdf5;
637 : #endif
638 :
639 3 : PetscFunctionBegin;
640 3 : PetscValidHeaderSpecific(nep,NEP_CLASSID,1);
641 3 : if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)nep),&viewer));
642 3 : PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
643 3 : PetscCheckSameComm(nep,1,viewer,2);
644 3 : NEPCheckSolved(nep,1);
645 3 : PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw));
646 3 : PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary));
647 : #if defined(PETSC_HAVE_HDF5)
648 : PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5));
649 : #endif
650 3 : PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
651 3 : if (isdraw) PetscCall(NEPValuesView_DRAW(nep,viewer));
652 3 : else if (isbinary) PetscCall(NEPValuesView_BINARY(nep,viewer));
653 : #if defined(PETSC_HAVE_HDF5)
654 : else if (ishdf5) PetscCall(NEPValuesView_HDF5(nep,viewer));
655 : #endif
656 2 : else if (isascii) {
657 2 : PetscCall(PetscViewerGetFormat(viewer,&format));
658 2 : switch (format) {
659 1 : case PETSC_VIEWER_DEFAULT:
660 : case PETSC_VIEWER_ASCII_INFO:
661 : case PETSC_VIEWER_ASCII_INFO_DETAIL:
662 1 : PetscCall(NEPValuesView_ASCII(nep,viewer));
663 : break;
664 1 : case PETSC_VIEWER_ASCII_MATLAB:
665 1 : PetscCall(NEPValuesView_MATLAB(nep,viewer));
666 : break;
667 0 : default:
668 0 : PetscCall(PetscInfo(nep,"Unsupported viewer format %s\n",PetscViewerFormats[format]));
669 : }
670 : }
671 3 : PetscFunctionReturn(PETSC_SUCCESS);
672 : }
673 :
674 : /*@
675 : NEPValuesViewFromOptions - Processes command line options to determine if/how
676 : the computed eigenvalues are to be viewed.
677 :
678 : Collective
679 :
680 : Input Parameter:
681 : . nep - the nonlinear eigensolver context
682 :
683 : Level: developer
684 :
685 : .seealso: NEPValuesView()
686 : @*/
687 160 : PetscErrorCode NEPValuesViewFromOptions(NEP nep)
688 : {
689 160 : PetscViewer viewer;
690 160 : PetscBool flg;
691 160 : static PetscBool incall = PETSC_FALSE;
692 160 : PetscViewerFormat format;
693 :
694 160 : PetscFunctionBegin;
695 160 : if (incall) PetscFunctionReturn(PETSC_SUCCESS);
696 160 : incall = PETSC_TRUE;
697 160 : PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)nep),((PetscObject)nep)->options,((PetscObject)nep)->prefix,"-nep_view_values",&viewer,&format,&flg));
698 160 : if (flg) {
699 3 : PetscCall(PetscViewerPushFormat(viewer,format));
700 3 : PetscCall(NEPValuesView(nep,viewer));
701 3 : PetscCall(PetscViewerPopFormat(viewer));
702 3 : PetscCall(PetscViewerDestroy(&viewer));
703 : }
704 160 : incall = PETSC_FALSE;
705 160 : PetscFunctionReturn(PETSC_SUCCESS);
706 : }
707 :
708 : /*@
709 : NEPVectorsView - Outputs computed eigenvectors to a viewer.
710 :
711 : Collective
712 :
713 : Input Parameters:
714 : + nep - the nonlinear eigensolver context
715 : - viewer - the viewer
716 :
717 : Options Database Key:
718 : . -nep_view_vectors - output eigenvectors.
719 :
720 : Notes:
721 : If PETSc was configured with real scalars, complex conjugate eigenvectors
722 : will be viewed as two separate real vectors, one containing the real part
723 : and another one containing the imaginary part.
724 :
725 : If left eigenvectors were computed with a two-sided eigensolver, the right
726 : and left eigenvectors are interleaved, that is, the vectors are output in
727 : the following order X0, Y0, X1, Y1, X2, Y2, ...
728 :
729 : Level: intermediate
730 :
731 : .seealso: NEPSolve(), NEPValuesView(), NEPErrorView()
732 : @*/
733 3 : PetscErrorCode NEPVectorsView(NEP nep,PetscViewer viewer)
734 : {
735 3 : PetscInt i,k;
736 3 : Vec xr,xi=NULL;
737 :
738 3 : PetscFunctionBegin;
739 3 : PetscValidHeaderSpecific(nep,NEP_CLASSID,1);
740 3 : if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)nep),&viewer));
741 3 : PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
742 3 : PetscCheckSameComm(nep,1,viewer,2);
743 3 : NEPCheckSolved(nep,1);
744 3 : if (nep->nconv) {
745 3 : PetscCall(NEPComputeVectors(nep));
746 3 : PetscCall(BVCreateVec(nep->V,&xr));
747 : #if !defined(PETSC_USE_COMPLEX)
748 : PetscCall(BVCreateVec(nep->V,&xi));
749 : #endif
750 6 : for (i=0;i<nep->nconv;i++) {
751 3 : k = nep->perm[i];
752 3 : PetscCall(BV_GetEigenvector(nep->V,k,nep->eigi[k],xr,xi));
753 3 : PetscCall(SlepcViewEigenvector(viewer,xr,xi,"X",i,(PetscObject)nep));
754 3 : if (nep->twosided) {
755 1 : PetscCall(BV_GetEigenvector(nep->W,k,nep->eigi[k],xr,xi));
756 3 : PetscCall(SlepcViewEigenvector(viewer,xr,xi,"Y",i,(PetscObject)nep));
757 : }
758 : }
759 3 : PetscCall(VecDestroy(&xr));
760 : #if !defined(PETSC_USE_COMPLEX)
761 : PetscCall(VecDestroy(&xi));
762 : #endif
763 : }
764 3 : PetscFunctionReturn(PETSC_SUCCESS);
765 : }
766 :
767 : /*@
768 : NEPVectorsViewFromOptions - Processes command line options to determine if/how
769 : the computed eigenvectors are to be viewed.
770 :
771 : Collective
772 :
773 : Input Parameter:
774 : . nep - the nonlinear eigensolver context
775 :
776 : Level: developer
777 :
778 : .seealso: NEPVectorsView()
779 : @*/
780 160 : PetscErrorCode NEPVectorsViewFromOptions(NEP nep)
781 : {
782 160 : PetscViewer viewer;
783 160 : PetscBool flg = PETSC_FALSE;
784 160 : static PetscBool incall = PETSC_FALSE;
785 160 : PetscViewerFormat format;
786 :
787 160 : PetscFunctionBegin;
788 160 : if (incall) PetscFunctionReturn(PETSC_SUCCESS);
789 160 : incall = PETSC_TRUE;
790 160 : PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)nep),((PetscObject)nep)->options,((PetscObject)nep)->prefix,"-nep_view_vectors",&viewer,&format,&flg));
791 160 : if (flg) {
792 3 : PetscCall(PetscViewerPushFormat(viewer,format));
793 3 : PetscCall(NEPVectorsView(nep,viewer));
794 3 : PetscCall(PetscViewerPopFormat(viewer));
795 3 : PetscCall(PetscViewerDestroy(&viewer));
796 : }
797 160 : incall = PETSC_FALSE;
798 160 : PetscFunctionReturn(PETSC_SUCCESS);
799 : }
|