Actual source code: nepmon.c
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
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 monitors
12: */
14: #include <slepc/private/nepimpl.h>
15: #include <petscdraw.h>
17: /*
18: Runs the user provided monitor routines, if any.
19: */
20: PetscErrorCode NEPMonitor(NEP nep,PetscInt it,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest)
21: {
22: PetscInt i,n = nep->numbermonitors;
24: PetscFunctionBegin;
25: for (i=0;i<n;i++) PetscCall((*nep->monitor[i])(nep,it,nconv,eigr,eigi,errest,nest,nep->monitorcontext[i]));
26: PetscFunctionReturn(PETSC_SUCCESS);
27: }
29: /*@C
30: NEPMonitorSet - Sets an ADDITIONAL function to be called at every
31: iteration to monitor the error estimates for each requested eigenpair.
33: Logically Collective
35: Input Parameters:
36: + nep - eigensolver context obtained from NEPCreate()
37: . monitor - pointer to function (if this is NULL, it turns off monitoring), see NEPMonitorFn
38: . mctx - [optional] context for private data for the
39: monitor routine (use NULL if no context is desired)
40: - monitordestroy - [optional] routine that frees monitor context (may be NULL),
41: see PetscCtxDestroyFn for the calling sequence
43: Options Database Keys:
44: + -nep_monitor - print only the first error estimate
45: . -nep_monitor_all - print error estimates at each iteration
46: . -nep_monitor_conv - print the eigenvalue approximations only when
47: convergence has been reached
48: . -nep_monitor draw::draw_lg - sets line graph monitor for the first unconverged
49: approximate eigenvalue
50: . -nep_monitor_all draw::draw_lg - sets line graph monitor for all unconverged
51: approximate eigenvalues
52: . -nep_monitor_conv draw::draw_lg - sets line graph monitor for convergence history
53: - -nep_monitor_cancel - cancels all monitors that have been hardwired into
54: a code by calls to NEPMonitorSet(), but does not cancel those set via
55: the options database.
57: Notes:
58: The options database option -nep_monitor and related options are the easiest way
59: to turn on NEP iteration monitoring.
61: NEPMonitorRegister() provides a way to associate an options database key with NEP
62: monitor function.
64: Several different monitoring routines may be set by calling NEPMonitorSet() multiple
65: times; all will be called in the order in which they were set.
67: Level: intermediate
69: .seealso: NEPMonitorFirst(), NEPMonitorAll(), NEPMonitorCancel()
70: @*/
71: PetscErrorCode NEPMonitorSet(NEP nep,NEPMonitorFn *monitor,void *mctx,PetscCtxDestroyFn *monitordestroy)
72: {
73: PetscInt i;
74: PetscBool identical;
76: PetscFunctionBegin;
78: for (i=0;i<nep->numbermonitors;i++) {
79: PetscCall(PetscMonitorCompare((PetscErrorCode(*)(void))(PetscVoidFn*)monitor,mctx,monitordestroy,(PetscErrorCode (*)(void))(PetscVoidFn*)nep->monitor[i],nep->monitorcontext[i],nep->monitordestroy[i],&identical));
80: if (identical) PetscFunctionReturn(PETSC_SUCCESS);
81: }
82: PetscCheck(nep->numbermonitors<MAXNEPMONITORS,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Too many NEP monitors set");
83: nep->monitor[nep->numbermonitors] = monitor;
84: nep->monitorcontext[nep->numbermonitors] = mctx;
85: nep->monitordestroy[nep->numbermonitors++] = monitordestroy;
86: PetscFunctionReturn(PETSC_SUCCESS);
87: }
89: /*@
90: NEPMonitorCancel - Clears all monitors for a NEP object.
92: Logically Collective
94: Input Parameters:
95: . nep - eigensolver context obtained from NEPCreate()
97: Options Database Key:
98: . -nep_monitor_cancel - Cancels all monitors that have been hardwired
99: into a code by calls to NEPMonitorSet(),
100: but does not cancel those set via the options database.
102: Level: intermediate
104: .seealso: NEPMonitorSet()
105: @*/
106: PetscErrorCode NEPMonitorCancel(NEP nep)
107: {
108: PetscInt i;
110: PetscFunctionBegin;
112: for (i=0; i<nep->numbermonitors; i++) {
113: if (nep->monitordestroy[i]) PetscCall((*nep->monitordestroy[i])(&nep->monitorcontext[i]));
114: }
115: nep->numbermonitors = 0;
116: PetscFunctionReturn(PETSC_SUCCESS);
117: }
119: /*@C
120: NEPGetMonitorContext - Gets the monitor context, as set by
121: NEPMonitorSet() for the FIRST monitor only.
123: Not Collective
125: Input Parameter:
126: . nep - eigensolver context obtained from NEPCreate()
128: Output Parameter:
129: . ctx - monitor context
131: Level: intermediate
133: .seealso: NEPMonitorSet(), NEPDefaultMonitor()
134: @*/
135: PetscErrorCode NEPGetMonitorContext(NEP nep,void *ctx)
136: {
137: PetscFunctionBegin;
139: *(void**)ctx = nep->monitorcontext[0];
140: PetscFunctionReturn(PETSC_SUCCESS);
141: }
143: /*@C
144: NEPMonitorFirst - Print the first unconverged approximate value and
145: error estimate at each iteration of the nonlinear eigensolver.
147: Collective
149: Input Parameters:
150: + nep - nonlinear eigensolver context
151: . its - iteration number
152: . nconv - number of converged eigenpairs so far
153: . eigr - real part of the eigenvalues
154: . eigi - imaginary part of the eigenvalues
155: . errest - error estimates
156: . nest - number of error estimates to display
157: - vf - viewer and format for monitoring
159: Options Database Key:
160: . -nep_monitor - activates NEPMonitorFirst()
162: Level: intermediate
164: .seealso: NEPMonitorSet(), NEPMonitorAll(), NEPMonitorConverged()
165: @*/
166: PetscErrorCode NEPMonitorFirst(NEP nep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
167: {
168: PetscViewer viewer = vf->viewer;
170: PetscFunctionBegin;
173: if (its==1 && ((PetscObject)nep)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer," Eigenvalue approximations and residual norms for %s solve.\n",((PetscObject)nep)->prefix));
174: if (nconv<nest) {
175: PetscCall(PetscViewerPushFormat(viewer,vf->format));
176: PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)nep)->tablevel));
177: PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " NEP nconv=%" PetscInt_FMT " first unconverged value (error)",its,nconv));
178: PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
179: #if defined(PETSC_USE_COMPLEX)
180: PetscCall(PetscViewerASCIIPrintf(viewer," %g%+gi",(double)PetscRealPart(eigr[nconv]),(double)PetscImaginaryPart(eigr[nconv])));
181: #else
182: PetscCall(PetscViewerASCIIPrintf(viewer," %g",(double)eigr[nconv]));
183: if (eigi[nconv]!=0.0) PetscCall(PetscViewerASCIIPrintf(viewer,"%+gi",(double)eigi[nconv]));
184: #endif
185: PetscCall(PetscViewerASCIIPrintf(viewer," (%10.8e)\n",(double)errest[nconv]));
186: PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
187: PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)nep)->tablevel));
188: PetscCall(PetscViewerPopFormat(viewer));
189: }
190: PetscFunctionReturn(PETSC_SUCCESS);
191: }
193: /*@C
194: NEPMonitorAll - Print the current approximate values and
195: error estimates at each iteration of the nonlinear eigensolver.
197: Collective
199: Input Parameters:
200: + nep - nonlinear eigensolver context
201: . its - iteration number
202: . nconv - number of converged eigenpairs so far
203: . eigr - real part of the eigenvalues
204: . eigi - imaginary part of the eigenvalues
205: . errest - error estimates
206: . nest - number of error estimates to display
207: - vf - viewer and format for monitoring
209: Options Database Key:
210: . -nep_monitor_all - activates NEPMonitorAll()
212: Level: intermediate
214: .seealso: NEPMonitorSet(), NEPMonitorFirst(), NEPMonitorConverged()
215: @*/
216: PetscErrorCode NEPMonitorAll(NEP nep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
217: {
218: PetscInt i;
219: PetscViewer viewer = vf->viewer;
221: PetscFunctionBegin;
224: PetscCall(PetscViewerPushFormat(viewer,vf->format));
225: PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)nep)->tablevel));
226: if (its==1 && ((PetscObject)nep)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer," Eigenvalue approximations and residual norms for %s solve.\n",((PetscObject)nep)->prefix));
227: PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " NEP nconv=%" PetscInt_FMT " Values (Errors)",its,nconv));
228: PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
229: for (i=0;i<nest;i++) {
230: #if defined(PETSC_USE_COMPLEX)
231: PetscCall(PetscViewerASCIIPrintf(viewer," %g%+gi",(double)PetscRealPart(eigr[i]),(double)PetscImaginaryPart(eigr[i])));
232: #else
233: PetscCall(PetscViewerASCIIPrintf(viewer," %g",(double)eigr[i]));
234: if (eigi[i]!=0.0) PetscCall(PetscViewerASCIIPrintf(viewer,"%+gi",(double)eigi[i]));
235: #endif
236: PetscCall(PetscViewerASCIIPrintf(viewer," (%10.8e)",(double)errest[i]));
237: }
238: PetscCall(PetscViewerASCIIPrintf(viewer,"\n"));
239: PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
240: PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)nep)->tablevel));
241: PetscCall(PetscViewerPopFormat(viewer));
242: PetscFunctionReturn(PETSC_SUCCESS);
243: }
245: /*@C
246: NEPMonitorConverged - Print the approximate values and
247: error estimates as they converge.
249: Collective
251: Input Parameters:
252: + nep - nonlinear eigensolver context
253: . its - iteration number
254: . nconv - number of converged eigenpairs so far
255: . eigr - real part of the eigenvalues
256: . eigi - imaginary part of the eigenvalues
257: . errest - error estimates
258: . nest - number of error estimates to display
259: - vf - viewer and format for monitoring
261: Options Database Key:
262: . -nep_monitor_conv - activates NEPMonitorConverged()
264: Level: intermediate
266: .seealso: NEPMonitorSet(), NEPMonitorFirst(), NEPMonitorAll()
267: @*/
268: PetscErrorCode NEPMonitorConverged(NEP nep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
269: {
270: PetscInt i;
271: PetscViewer viewer = vf->viewer;
272: SlepcConvMon ctx;
274: PetscFunctionBegin;
277: ctx = (SlepcConvMon)vf->data;
278: if (its==1 && ((PetscObject)nep)->prefix) PetscCall(PetscViewerASCIIPrintf(viewer," Convergence history for %s solve.\n",((PetscObject)nep)->prefix));
279: if (its==1) ctx->oldnconv = 0;
280: if (ctx->oldnconv!=nconv) {
281: PetscCall(PetscViewerPushFormat(viewer,vf->format));
282: PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)nep)->tablevel));
283: for (i=ctx->oldnconv;i<nconv;i++) {
284: PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " NEP converged value (error) #%" PetscInt_FMT,its,i));
285: PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
286: #if defined(PETSC_USE_COMPLEX)
287: PetscCall(PetscViewerASCIIPrintf(viewer," %g%+gi",(double)PetscRealPart(eigr[i]),(double)PetscImaginaryPart(eigr[i])));
288: #else
289: PetscCall(PetscViewerASCIIPrintf(viewer," %g",(double)eigr[i]));
290: if (eigi[i]!=0.0) PetscCall(PetscViewerASCIIPrintf(viewer,"%+gi",(double)eigi[i]));
291: #endif
292: PetscCall(PetscViewerASCIIPrintf(viewer," (%10.8e)\n",(double)errest[i]));
293: PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
294: }
295: PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)nep)->tablevel));
296: PetscCall(PetscViewerPopFormat(viewer));
297: ctx->oldnconv = nconv;
298: }
299: PetscFunctionReturn(PETSC_SUCCESS);
300: }
302: PetscErrorCode NEPMonitorConvergedCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
303: {
304: SlepcConvMon mctx;
306: PetscFunctionBegin;
307: PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
308: PetscCall(PetscNew(&mctx));
309: mctx->ctx = ctx;
310: (*vf)->data = (void*)mctx;
311: PetscFunctionReturn(PETSC_SUCCESS);
312: }
314: PetscErrorCode NEPMonitorConvergedDestroy(PetscViewerAndFormat **vf)
315: {
316: PetscFunctionBegin;
317: if (!*vf) PetscFunctionReturn(PETSC_SUCCESS);
318: PetscCall(PetscFree((*vf)->data));
319: PetscCall(PetscViewerDestroy(&(*vf)->viewer));
320: PetscCall(PetscFree(*vf));
321: PetscFunctionReturn(PETSC_SUCCESS);
322: }
324: /*@C
325: NEPMonitorFirstDrawLG - Plots the error estimate of the first unconverged
326: approximation at each iteration of the nonlinear eigensolver.
328: Collective
330: Input Parameters:
331: + nep - nonlinear eigensolver context
332: . its - iteration number
333: . nconv - number of converged eigenpairs so far
334: . eigr - real part of the eigenvalues
335: . eigi - imaginary part of the eigenvalues
336: . errest - error estimates
337: . nest - number of error estimates to display
338: - vf - viewer and format for monitoring
340: Options Database Key:
341: . -nep_monitor draw::draw_lg - activates NEPMonitorFirstDrawLG()
343: Level: intermediate
345: .seealso: NEPMonitorSet()
346: @*/
347: PetscErrorCode NEPMonitorFirstDrawLG(NEP nep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
348: {
349: PetscViewer viewer = vf->viewer;
350: PetscDrawLG lg;
351: PetscReal x,y;
353: PetscFunctionBegin;
356: PetscCall(PetscViewerPushFormat(viewer,vf->format));
357: PetscCall(PetscViewerDrawGetDrawLG(viewer,0,&lg));
358: if (its==1) {
359: PetscCall(PetscDrawLGReset(lg));
360: PetscCall(PetscDrawLGSetDimension(lg,1));
361: PetscCall(PetscDrawLGSetLimits(lg,1,1,PetscLog10Real(nep->tol)-2,0.0));
362: }
363: if (nconv<nest) {
364: x = (PetscReal)its;
365: if (errest[nconv] > 0.0) y = PetscLog10Real(errest[nconv]);
366: else y = 0.0;
367: PetscCall(PetscDrawLGAddPoint(lg,&x,&y));
368: if (its <= 20 || !(its % 5) || nep->reason) {
369: PetscCall(PetscDrawLGDraw(lg));
370: PetscCall(PetscDrawLGSave(lg));
371: }
372: }
373: PetscCall(PetscViewerPopFormat(viewer));
374: PetscFunctionReturn(PETSC_SUCCESS);
375: }
377: /*@C
378: NEPMonitorFirstDrawLGCreate - Creates the plotter for the first error estimate.
380: Collective
382: Input Parameters:
383: + viewer - the viewer
384: . format - the viewer format
385: - ctx - an optional user context
387: Output Parameter:
388: . vf - the viewer and format context
390: Level: intermediate
392: .seealso: NEPMonitorSet()
393: @*/
394: PetscErrorCode NEPMonitorFirstDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
395: {
396: PetscFunctionBegin;
397: PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
398: (*vf)->data = ctx;
399: PetscCall(PetscViewerMonitorLGSetUp(viewer,NULL,"First Error Estimate","Log Error Estimate",1,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300));
400: PetscFunctionReturn(PETSC_SUCCESS);
401: }
403: /*@C
404: NEPMonitorAllDrawLG - Plots the error estimate of all unconverged
405: approximations at each iteration of the nonlinear eigensolver.
407: Collective
409: Input Parameters:
410: + nep - nonlinear eigensolver context
411: . its - iteration number
412: . nconv - number of converged eigenpairs so far
413: . eigr - real part of the eigenvalues
414: . eigi - imaginary part of the eigenvalues
415: . errest - error estimates
416: . nest - number of error estimates to display
417: - vf - viewer and format for monitoring
419: Options Database Key:
420: . -nep_monitor_all draw::draw_lg - activates NEPMonitorAllDrawLG()
422: Level: intermediate
424: .seealso: NEPMonitorSet()
425: @*/
426: PetscErrorCode NEPMonitorAllDrawLG(NEP nep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
427: {
428: PetscViewer viewer = vf->viewer;
429: PetscDrawLG lg;
430: PetscInt i,n = PetscMin(nep->nev,255);
431: PetscReal *x,*y;
433: PetscFunctionBegin;
436: PetscCall(PetscViewerPushFormat(viewer,vf->format));
437: PetscCall(PetscViewerDrawGetDrawLG(viewer,0,&lg));
438: if (its==1) {
439: PetscCall(PetscDrawLGReset(lg));
440: PetscCall(PetscDrawLGSetDimension(lg,n));
441: PetscCall(PetscDrawLGSetLimits(lg,1,1,PetscLog10Real(nep->tol)-2,0.0));
442: }
443: PetscCall(PetscMalloc2(n,&x,n,&y));
444: for (i=0;i<n;i++) {
445: x[i] = (PetscReal)its;
446: if (i < nest && errest[i] > 0.0) y[i] = PetscLog10Real(errest[i]);
447: else y[i] = 0.0;
448: }
449: PetscCall(PetscDrawLGAddPoint(lg,x,y));
450: if (its <= 20 || !(its % 5) || nep->reason) {
451: PetscCall(PetscDrawLGDraw(lg));
452: PetscCall(PetscDrawLGSave(lg));
453: }
454: PetscCall(PetscFree2(x,y));
455: PetscCall(PetscViewerPopFormat(viewer));
456: PetscFunctionReturn(PETSC_SUCCESS);
457: }
459: /*@C
460: NEPMonitorAllDrawLGCreate - Creates the plotter for all the error estimates.
462: Collective
464: Input Parameters:
465: + viewer - the viewer
466: . format - the viewer format
467: - ctx - an optional user context
469: Output Parameter:
470: . vf - the viewer and format context
472: Level: intermediate
474: .seealso: NEPMonitorSet()
475: @*/
476: PetscErrorCode NEPMonitorAllDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
477: {
478: PetscFunctionBegin;
479: PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
480: (*vf)->data = ctx;
481: PetscCall(PetscViewerMonitorLGSetUp(viewer,NULL,"All Error Estimates","Log Error Estimate",1,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300));
482: PetscFunctionReturn(PETSC_SUCCESS);
483: }
485: /*@C
486: NEPMonitorConvergedDrawLG - Plots the number of converged eigenvalues
487: at each iteration of the nonlinear eigensolver.
489: Collective
491: Input Parameters:
492: + nep - nonlinear eigensolver context
493: . its - iteration number
494: . nconv - number of converged eigenpairs so far
495: . eigr - real part of the eigenvalues
496: . eigi - imaginary part of the eigenvalues
497: . errest - error estimates
498: . nest - number of error estimates to display
499: - vf - viewer and format for monitoring
501: Options Database Key:
502: . -nep_monitor_conv draw::draw_lg - activates NEPMonitorConvergedDrawLG()
504: Level: intermediate
506: .seealso: NEPMonitorSet()
507: @*/
508: PetscErrorCode NEPMonitorConvergedDrawLG(NEP nep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *vf)
509: {
510: PetscViewer viewer = vf->viewer;
511: PetscDrawLG lg;
512: PetscReal x,y;
514: PetscFunctionBegin;
517: PetscCall(PetscViewerPushFormat(viewer,vf->format));
518: PetscCall(PetscViewerDrawGetDrawLG(viewer,0,&lg));
519: if (its==1) {
520: PetscCall(PetscDrawLGReset(lg));
521: PetscCall(PetscDrawLGSetDimension(lg,1));
522: PetscCall(PetscDrawLGSetLimits(lg,1,1,0,nep->nev));
523: }
524: x = (PetscReal)its;
525: y = (PetscReal)nep->nconv;
526: PetscCall(PetscDrawLGAddPoint(lg,&x,&y));
527: if (its <= 20 || !(its % 5) || nep->reason) {
528: PetscCall(PetscDrawLGDraw(lg));
529: PetscCall(PetscDrawLGSave(lg));
530: }
531: PetscCall(PetscViewerPopFormat(viewer));
532: PetscFunctionReturn(PETSC_SUCCESS);
533: }
535: /*@C
536: NEPMonitorConvergedDrawLGCreate - Creates the plotter for the convergence history.
538: Collective
540: Input Parameters:
541: + viewer - the viewer
542: . format - the viewer format
543: - ctx - an optional user context
545: Output Parameter:
546: . vf - the viewer and format context
548: Level: intermediate
550: .seealso: NEPMonitorSet()
551: @*/
552: PetscErrorCode NEPMonitorConvergedDrawLGCreate(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **vf)
553: {
554: SlepcConvMon mctx;
556: PetscFunctionBegin;
557: PetscCall(PetscViewerAndFormatCreate(viewer,format,vf));
558: PetscCall(PetscNew(&mctx));
559: mctx->ctx = ctx;
560: (*vf)->data = (void*)mctx;
561: PetscCall(PetscViewerMonitorLGSetUp(viewer,NULL,"Convergence History","Number Converged",1,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300));
562: PetscFunctionReturn(PETSC_SUCCESS);
563: }