Actual source code: svdbasic.c
1: /*
2: The basic SVD routines, Create, View, etc. are here.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2011, Universitat Politecnica de Valencia, Spain
8: This file is part of SLEPc.
9:
10: SLEPc is free software: you can redistribute it and/or modify it under the
11: terms of version 3 of the GNU Lesser General Public License as published by
12: the Free Software Foundation.
14: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
15: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
16: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
17: more details.
19: You should have received a copy of the GNU Lesser General Public License
20: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
21: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
22: */
24: #include <private/svdimpl.h> /*I "slepcsvd.h" I*/
26: PetscFList SVDList = 0;
27: PetscBool SVDRegisterAllCalled = PETSC_FALSE;
28: PetscClassId SVD_CLASSID = 0;
29: PetscLogEvent SVD_SetUp = 0,SVD_Solve = 0,SVD_Dense = 0;
30: static PetscBool SVDPackageInitialized = PETSC_FALSE;
34: /*@C
35: SVDFinalizePackage - This function destroys everything in the Slepc interface
36: to the SVD package. It is called from SlepcFinalize().
38: Level: developer
40: .seealso: SlepcFinalize()
41: @*/
42: PetscErrorCode SVDFinalizePackage(void)
43: {
45: SVDPackageInitialized = PETSC_FALSE;
46: SVDList = 0;
47: SVDRegisterAllCalled = PETSC_FALSE;
48: return(0);
49: }
53: /*@C
54: SVDInitializePackage - This function initializes everything in the SVD package. It is called
55: from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to SVDCreate()
56: when using static libraries.
58: Input Parameter:
59: . path - The dynamic library path, or PETSC_NULL
61: Level: developer
63: .seealso: SlepcInitialize()
64: @*/
65: PetscErrorCode SVDInitializePackage(const char *path)
66: {
67: char logList[256];
68: char *className;
69: PetscBool opt;
73: if (SVDPackageInitialized) return(0);
74: SVDPackageInitialized = PETSC_TRUE;
75: /* Register Classes */
76: PetscClassIdRegister("Singular Value Solver",&SVD_CLASSID);
77: /* Register Constructors */
78: SVDRegisterAll(path);
79: /* Register Events */
80: PetscLogEventRegister("SVDSetUp",SVD_CLASSID,&SVD_SetUp);
81: PetscLogEventRegister("SVDSolve",SVD_CLASSID,&SVD_Solve);
82: PetscLogEventRegister("SVDDense",SVD_CLASSID,&SVD_Dense);
83: /* Process info exclusions */
84: PetscOptionsGetString(PETSC_NULL,"-info_exclude",logList,256,&opt);
85: if (opt) {
86: PetscStrstr(logList,"svd",&className);
87: if (className) {
88: PetscInfoDeactivateClass(SVD_CLASSID);
89: }
90: }
91: /* Process summary exclusions */
92: PetscOptionsGetString(PETSC_NULL,"-log_summary_exclude",logList,256,&opt);
93: if (opt) {
94: PetscStrstr(logList,"svd",&className);
95: if (className) {
96: PetscLogEventDeactivateClass(SVD_CLASSID);
97: }
98: }
99: PetscRegisterFinalize(SVDFinalizePackage);
100: return(0);
101: }
105: /*@C
106: SVDView - Prints the SVD data structure.
108: Collective on SVD
110: Input Parameters:
111: + svd - the singular value solver context
112: - viewer - optional visualization context
114: Options Database Key:
115: . -svd_view - Calls SVDView() at end of SVDSolve()
117: Note:
118: The available visualization contexts include
119: + PETSC_VIEWER_STDOUT_SELF - standard output (default)
120: - PETSC_VIEWER_STDOUT_WORLD - synchronized standard
121: output where only the first processor opens
122: the file. All other processors send their
123: data to the first processor to print.
125: The user can open an alternative visualization context with
126: PetscViewerASCIIOpen() - output to a specified file.
128: Level: beginner
130: .seealso: STView(), PetscViewerASCIIOpen()
131: @*/
132: PetscErrorCode SVDView(SVD svd,PetscViewer viewer)
133: {
135: PetscBool isascii;
139: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(((PetscObject)svd)->comm);
143: PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
144: if (isascii) {
145: PetscObjectPrintClassNamePrefixType((PetscObject)svd,viewer,"SVD Object");
146: if (svd->ops->view) {
147: PetscViewerASCIIPushTab(viewer);
148: (*svd->ops->view)(svd,viewer);
149: PetscViewerASCIIPopTab(viewer);
150: }
151: switch (svd->transmode) {
152: case SVD_TRANSPOSE_EXPLICIT:
153: PetscViewerASCIIPrintf(viewer," transpose mode: explicit\n");
154: break;
155: case SVD_TRANSPOSE_IMPLICIT:
156: PetscViewerASCIIPrintf(viewer," transpose mode: implicit\n");
157: break;
158: default:
159: PetscViewerASCIIPrintf(viewer," transpose mode: not yet set\n");
160: }
161: if (svd->which == SVD_LARGEST) {
162: PetscViewerASCIIPrintf(viewer," selected portion of the spectrum: largest\n");
163: } else {
164: PetscViewerASCIIPrintf(viewer," selected portion of the spectrum: smallest\n");
165: }
166: PetscViewerASCIIPrintf(viewer," number of singular values (nsv): %D\n",svd->nsv);
167: PetscViewerASCIIPrintf(viewer," number of column vectors (ncv): %D\n",svd->ncv);
168: PetscViewerASCIIPrintf(viewer," maximum dimension of projected problem (mpd): %D\n",svd->mpd);
169: PetscViewerASCIIPrintf(viewer," maximum number of iterations: %D\n",svd->max_it);
170: PetscViewerASCIIPrintf(viewer," tolerance: %G\n",svd->tol);
171: if (svd->nini!=0) {
172: PetscViewerASCIIPrintf(viewer," dimension of user-provided initial space: %D\n",PetscAbs(svd->nini));
173: }
174: } else {
175: if (svd->ops->view) {
176: (*svd->ops->view)(svd,viewer);
177: }
178: }
179: if (!svd->ip) { SVDGetIP(svd,&svd->ip); }
180: IPView(svd->ip,viewer);
181: return(0);
182: }
186: /*@
187: SVDPrintSolution - Prints the computed singular values.
189: Collective on SVD
191: Input Parameters:
192: + svd - the singular value solver context
193: - viewer - optional visualization context
195: Options Database:
196: . -svd_terse - print only minimal information
198: Note:
199: By default, this function prints a table with singular values and associated
200: relative errors. With -svd_terse only the singular values are printed.
202: Level: intermediate
204: .seealso: PetscViewerASCIIOpen()
205: @*/
206: PetscErrorCode SVDPrintSolution(SVD svd,PetscViewer viewer)
207: {
208: PetscBool terse,errok,isascii;
209: PetscReal error,sigma;
210: PetscInt i,j;
215: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(((PetscObject)svd)->comm);
218: if (!svd->sigma) {
219: SETERRQ(((PetscObject)svd)->comm,PETSC_ERR_ARG_WRONGSTATE,"SVDSolve must be called first");
220: }
221: PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
222: if (!isascii) return(0);
224: PetscOptionsHasName(PETSC_NULL,"-svd_terse",&terse);
225: if (terse) {
226: if (svd->nconv<svd->nsv) {
227: PetscViewerASCIIPrintf(viewer," Problem: less than %D singular values converged\n\n",svd->nsv);
228: } else {
229: errok = PETSC_TRUE;
230: for (i=0;i<svd->nsv;i++) {
231: SVDComputeRelativeError(svd,i,&error);
232: errok = (errok && error<svd->tol)? PETSC_TRUE: PETSC_FALSE;
233: }
234: if (errok) {
235: PetscViewerASCIIPrintf(viewer," All requested singular values computed up to the required tolerance:");
236: for (i=0;i<=(svd->nsv-1)/8;i++) {
237: PetscViewerASCIIPrintf(viewer,"\n ");
238: for (j=0;j<PetscMin(8,svd->nsv-8*i);j++) {
239: SVDGetSingularTriplet(svd,8*i+j,&sigma,PETSC_NULL,PETSC_NULL);
240: PetscViewerASCIIPrintf(viewer,"%.5F",sigma);
241: if (8*i+j+1<svd->nsv) { PetscViewerASCIIPrintf(viewer,", "); }
242: }
243: }
244: PetscViewerASCIIPrintf(viewer,"\n\n");
245: } else {
246: PetscViewerASCIIPrintf(viewer," Problem: some of the first %D relative errors are higher than the tolerance\n\n",svd->nsv);
247: }
248: }
249: } else {
250: PetscViewerASCIIPrintf(viewer," Number of converged approximate singular triplets: %D\n\n",svd->nconv);
251: if (svd->nconv>0) {
252: PetscViewerASCIIPrintf(viewer,
253: " sigma relative error\n"
254: " --------------------- ------------------\n");
255: for (i=0;i<svd->nconv;i++) {
256: SVDGetSingularTriplet(svd,i,&sigma,PETSC_NULL,PETSC_NULL);
257: SVDComputeRelativeError(svd,i,&error);
258: PetscViewerASCIIPrintf(viewer," % 6F %12G\n",sigma,error);
259: }
260: PetscViewerASCIIPrintf(viewer,"\n");
261: }
262: }
263: return(0);
264: }
268: /*@C
269: SVDCreate - Creates the default SVD context.
271: Collective on MPI_Comm
273: Input Parameter:
274: . comm - MPI communicator
276: Output Parameter:
277: . svd - location to put the SVD context
279: Note:
280: The default SVD type is SVDCROSS
282: Level: beginner
284: .seealso: SVDSetUp(), SVDSolve(), SVDDestroy(), SVD
285: @*/
286: PetscErrorCode SVDCreate(MPI_Comm comm,SVD *outsvd)
287: {
289: SVD svd;
293: *outsvd = 0;
294: PetscHeaderCreate(svd,_p_SVD,struct _SVDOps,SVD_CLASSID,-1,"SVD","Singular Value Decomposition","SVD",comm,SVDDestroy,SVDView);
296: svd->OP = PETSC_NULL;
297: svd->A = PETSC_NULL;
298: svd->AT = PETSC_NULL;
299: svd->transmode = (SVDTransposeMode)PETSC_DECIDE;
300: svd->sigma = PETSC_NULL;
301: svd->perm = PETSC_NULL;
302: svd->U = PETSC_NULL;
303: svd->V = PETSC_NULL;
304: svd->IS = PETSC_NULL;
305: svd->tl = PETSC_NULL;
306: svd->tr = PETSC_NULL;
307: svd->rand = PETSC_NULL;
308: svd->which = SVD_LARGEST;
309: svd->n = 0;
310: svd->nconv = 0;
311: svd->nsv = 1;
312: svd->ncv = 0;
313: svd->mpd = 0;
314: svd->nini = 0;
315: svd->its = 0;
316: svd->max_it = 0;
317: svd->tol = PETSC_DEFAULT;
318: svd->errest = PETSC_NULL;
319: svd->data = PETSC_NULL;
320: svd->setupcalled = 0;
321: svd->reason = SVD_CONVERGED_ITERATING;
322: svd->numbermonitors = 0;
323: svd->matvecs = 0;
324: svd->trackall = PETSC_FALSE;
326: PetscRandomCreate(comm,&svd->rand);
327: PetscLogObjectParent(svd,svd->rand);
328: *outsvd = svd;
329: return(0);
330: }
331:
334: /*@
335: SVDReset - Resets the SVD context to the setupcalled=0 state and removes any
336: allocated objects.
338: Collective on SVD
340: Input Parameter:
341: . svd - singular value solver context obtained from SVDCreate()
343: Level: advanced
345: .seealso: SVDDestroy()
346: @*/
347: PetscErrorCode SVDReset(SVD svd)
348: {
353: if (svd->ops->reset) { (svd->ops->reset)(svd); }
354: if (svd->ip) { IPReset(svd->ip); }
355: MatDestroy(&svd->OP);
356: MatDestroy(&svd->A);
357: MatDestroy(&svd->AT);
358: VecDestroy(&svd->tl);
359: VecDestroy(&svd->tr);
360: if (svd->n) {
361: PetscFree(svd->sigma);
362: PetscFree(svd->perm);
363: PetscFree(svd->errest);
364: VecDestroyVecs(svd->n,&svd->U);
365: VecDestroyVecs(svd->n,&svd->V);
366: }
367: svd->transmode = (SVDTransposeMode)PETSC_DECIDE;
368: svd->matvecs = 0;
369: svd->setupcalled = 0;
370: return(0);
371: }
375: /*@C
376: SVDDestroy - Destroys the SVD context.
378: Collective on SVD
380: Input Parameter:
381: . svd - singular value solver context obtained from SVDCreate()
383: Level: beginner
385: .seealso: SVDCreate(), SVDSetUp(), SVDSolve()
386: @*/
387: PetscErrorCode SVDDestroy(SVD *svd)
388: {
390:
392: if (!*svd) return(0);
394: if (--((PetscObject)(*svd))->refct > 0) { *svd = 0; return(0); }
395: SVDReset(*svd);
396: PetscObjectDepublish(*svd);
397: if ((*svd)->ops->destroy) { (*(*svd)->ops->destroy)(*svd); }
398: IPDestroy(&(*svd)->ip);
399: PetscRandomDestroy(&(*svd)->rand);
400: SVDMonitorCancel(*svd);
401: PetscHeaderDestroy(svd);
402: return(0);
403: }
407: /*@C
408: SVDSetType - Selects the particular solver to be used in the SVD object.
410: Logically Collective on SVD
412: Input Parameters:
413: + svd - the singular value solver context
414: - type - a known method
416: Options Database Key:
417: . -svd_type <method> - Sets the method; use -help for a list
418: of available methods
419:
420: Notes:
421: See "slepc/include/slepcsvd.h" for available methods. The default
422: is SVDCROSS.
424: Normally, it is best to use the SVDSetFromOptions() command and
425: then set the SVD type from the options database rather than by using
426: this routine. Using the options database provides the user with
427: maximum flexibility in evaluating the different available methods.
428: The SVDSetType() routine is provided for those situations where it
429: is necessary to set the iterative solver independently of the command
430: line or options database.
432: Level: intermediate
434: .seealso: SVDType
435: @*/
436: PetscErrorCode SVDSetType(SVD svd,const SVDType type)
437: {
438: PetscErrorCode ierr,(*r)(SVD);
439: PetscBool match;
445: PetscTypeCompare((PetscObject)svd,type,&match);
446: if (match) return(0);
448: PetscFListFind(SVDList,((PetscObject)svd)->comm,type,PETSC_TRUE,(void (**)(void)) &r);
449: if (!r) SETERRQ1(((PetscObject)svd)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown SVD type given: %s",type);
451: if (svd->ops->destroy) { (*svd->ops->destroy)(svd); }
452: PetscMemzero(svd->ops,sizeof(struct _SVDOps));
454: svd->setupcalled = 0;
455: PetscObjectChangeTypeName((PetscObject)svd,type);
456: (*r)(svd);
457: return(0);
458: }
462: /*@C
463: SVDGetType - Gets the SVD type as a string from the SVD object.
465: Not Collective
467: Input Parameter:
468: . svd - the singular value solver context
470: Output Parameter:
471: . name - name of SVD method
473: Level: intermediate
475: .seealso: SVDSetType()
476: @*/
477: PetscErrorCode SVDGetType(SVD svd,const SVDType *type)
478: {
482: *type = ((PetscObject)svd)->type_name;
483: return(0);
484: }
488: /*@C
489: SVDRegister - See SVDRegisterDynamic()
491: Level: advanced
492: @*/
493: PetscErrorCode SVDRegister(const char *sname,const char *path,const char *name,PetscErrorCode (*function)(SVD))
494: {
496: char fullname[PETSC_MAX_PATH_LEN];
499: PetscFListConcat(path,name,fullname);
500: PetscFListAdd(&SVDList,sname,fullname,(void (*)(void))function);
501: return(0);
502: }
506: /*@
507: SVDRegisterDestroy - Frees the list of SVD methods that were
508: registered by SVDRegisterDynamic().
510: Not Collective
512: Level: advanced
514: .seealso: SVDRegisterDynamic(), SVDRegisterAll()
515: @*/
516: PetscErrorCode SVDRegisterDestroy(void)
517: {
521: PetscFListDestroy(&SVDList);
522: SVDRegisterAllCalled = PETSC_FALSE;
523: return(0);
524: }
528: /*@
529: SVDSetIP - Associates an inner product object to the
530: singular value solver.
532: Collective on SVD
534: Input Parameters:
535: + svd - singular value solver context obtained from SVDCreate()
536: - ip - the inner product object
538: Note:
539: Use SVDGetIP() to retrieve the inner product context (for example,
540: to free it at the end of the computations).
542: Level: advanced
544: .seealso: SVDGetIP()
545: @*/
546: PetscErrorCode SVDSetIP(SVD svd,IP ip)
547: {
554: PetscObjectReference((PetscObject)ip);
555: IPDestroy(&svd->ip);
556: svd->ip = ip;
557: PetscLogObjectParent(svd,svd->ip);
558: return(0);
559: }
563: /*@C
564: SVDGetIP - Obtain the inner product object associated
565: to the singular value solver object.
567: Not Collective
569: Input Parameters:
570: . svd - singular value solver context obtained from SVDCreate()
572: Output Parameter:
573: . ip - inner product context
575: Level: advanced
577: .seealso: SVDSetIP()
578: @*/
579: PetscErrorCode SVDGetIP(SVD svd,IP *ip)
580: {
586: if (!svd->ip) {
587: IPCreate(((PetscObject)svd)->comm,&svd->ip);
588: PetscLogObjectParent(svd,svd->ip);
589: }
590: *ip = svd->ip;
591: return(0);
592: }