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-2012, 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 <slepc-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;
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: /* Process info exclusions */
83: PetscOptionsGetString(PETSC_NULL,"-info_exclude",logList,256,&opt);
84: if (opt) {
85: PetscStrstr(logList,"svd",&className);
86: if (className) {
87: PetscInfoDeactivateClass(SVD_CLASSID);
88: }
89: }
90: /* Process summary exclusions */
91: PetscOptionsGetString(PETSC_NULL,"-log_summary_exclude",logList,256,&opt);
92: if (opt) {
93: PetscStrstr(logList,"svd",&className);
94: if (className) {
95: PetscLogEventDeactivateClass(SVD_CLASSID);
96: }
97: }
98: PetscRegisterFinalize(SVDFinalizePackage);
99: return(0);
100: }
104: /*@C
105: SVDView - Prints the SVD data structure.
107: Collective on SVD
109: Input Parameters:
110: + svd - the singular value solver context
111: - viewer - optional visualization context
113: Options Database Key:
114: . -svd_view - Calls SVDView() at end of SVDSolve()
116: Note:
117: The available visualization contexts include
118: + PETSC_VIEWER_STDOUT_SELF - standard output (default)
119: - PETSC_VIEWER_STDOUT_WORLD - synchronized standard
120: output where only the first processor opens
121: the file. All other processors send their
122: data to the first processor to print.
124: The user can open an alternative visualization context with
125: PetscViewerASCIIOpen() - output to a specified file.
127: Level: beginner
129: .seealso: STView(), PetscViewerASCIIOpen()
130: @*/
131: PetscErrorCode SVDView(SVD svd,PetscViewer viewer)
132: {
134: PetscBool isascii,isshell;
138: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(((PetscObject)svd)->comm);
142: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
143: if (isascii) {
144: PetscObjectPrintClassNamePrefixType((PetscObject)svd,viewer,"SVD Object");
145: if (svd->ops->view) {
146: PetscViewerASCIIPushTab(viewer);
147: (*svd->ops->view)(svd,viewer);
148: PetscViewerASCIIPopTab(viewer);
149: }
150: switch (svd->transmode) {
151: case SVD_TRANSPOSE_EXPLICIT:
152: PetscViewerASCIIPrintf(viewer," transpose mode: explicit\n");
153: break;
154: case SVD_TRANSPOSE_IMPLICIT:
155: PetscViewerASCIIPrintf(viewer," transpose mode: implicit\n");
156: break;
157: default:
158: PetscViewerASCIIPrintf(viewer," transpose mode: not yet set\n");
159: }
160: if (svd->which == SVD_LARGEST) {
161: PetscViewerASCIIPrintf(viewer," selected portion of the spectrum: largest\n");
162: } else {
163: PetscViewerASCIIPrintf(viewer," selected portion of the spectrum: smallest\n");
164: }
165: PetscViewerASCIIPrintf(viewer," number of singular values (nsv): %D\n",svd->nsv);
166: PetscViewerASCIIPrintf(viewer," number of column vectors (ncv): %D\n",svd->ncv);
167: PetscViewerASCIIPrintf(viewer," maximum dimension of projected problem (mpd): %D\n",svd->mpd);
168: PetscViewerASCIIPrintf(viewer," maximum number of iterations: %D\n",svd->max_it);
169: PetscViewerASCIIPrintf(viewer," tolerance: %G\n",svd->tol);
170: if (svd->nini!=0) {
171: PetscViewerASCIIPrintf(viewer," dimension of user-provided initial space: %D\n",PetscAbs(svd->nini));
172: }
173: } else {
174: if (svd->ops->view) {
175: (*svd->ops->view)(svd,viewer);
176: }
177: }
178: PetscObjectTypeCompareAny((PetscObject)svd,&isshell,SVDCROSS,SVDCYCLIC,"");
179: if (!isshell) {
180: if (!svd->ip) { SVDGetIP(svd,&svd->ip); }
181: IPView(svd->ip,viewer);
182: if (!svd->ds) { SVDGetDS(svd,&svd->ds); }
183: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
184: DSView(svd->ds,viewer);
185: PetscViewerPopFormat(viewer);
186: }
187: return(0);
188: }
192: /*@
193: SVDPrintSolution - Prints the computed singular values.
195: Collective on SVD
197: Input Parameters:
198: + svd - the singular value solver context
199: - viewer - optional visualization context
201: Options Database:
202: . -svd_terse - print only minimal information
204: Note:
205: By default, this function prints a table with singular values and associated
206: relative errors. With -svd_terse only the singular values are printed.
208: Level: intermediate
210: .seealso: PetscViewerASCIIOpen()
211: @*/
212: PetscErrorCode SVDPrintSolution(SVD svd,PetscViewer viewer)
213: {
214: PetscBool terse,errok,isascii;
215: PetscReal error,sigma;
216: PetscInt i,j;
221: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(((PetscObject)svd)->comm);
224: if (!svd->sigma) SETERRQ(((PetscObject)svd)->comm,PETSC_ERR_ARG_WRONGSTATE,"SVDSolve must be called first");
225: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
226: if (!isascii) return(0);
228: PetscOptionsHasName(PETSC_NULL,"-svd_terse",&terse);
229: if (terse) {
230: if (svd->nconv<svd->nsv) {
231: PetscViewerASCIIPrintf(viewer," Problem: less than %D singular values converged\n\n",svd->nsv);
232: } else {
233: errok = PETSC_TRUE;
234: for (i=0;i<svd->nsv;i++) {
235: SVDComputeRelativeError(svd,i,&error);
236: errok = (errok && error<svd->tol)? PETSC_TRUE: PETSC_FALSE;
237: }
238: if (errok) {
239: PetscViewerASCIIPrintf(viewer," All requested singular values computed up to the required tolerance:");
240: for (i=0;i<=(svd->nsv-1)/8;i++) {
241: PetscViewerASCIIPrintf(viewer,"\n ");
242: for (j=0;j<PetscMin(8,svd->nsv-8*i);j++) {
243: SVDGetSingularTriplet(svd,8*i+j,&sigma,PETSC_NULL,PETSC_NULL);
244: PetscViewerASCIIPrintf(viewer,"%.5F",sigma);
245: if (8*i+j+1<svd->nsv) { PetscViewerASCIIPrintf(viewer,", "); }
246: }
247: }
248: PetscViewerASCIIPrintf(viewer,"\n\n");
249: } else {
250: PetscViewerASCIIPrintf(viewer," Problem: some of the first %D relative errors are higher than the tolerance\n\n",svd->nsv);
251: }
252: }
253: } else {
254: PetscViewerASCIIPrintf(viewer," Number of converged approximate singular triplets: %D\n\n",svd->nconv);
255: if (svd->nconv>0) {
256: PetscViewerASCIIPrintf(viewer,
257: " sigma relative error\n"
258: " --------------------- ------------------\n");
259: for (i=0;i<svd->nconv;i++) {
260: SVDGetSingularTriplet(svd,i,&sigma,PETSC_NULL,PETSC_NULL);
261: SVDComputeRelativeError(svd,i,&error);
262: PetscViewerASCIIPrintf(viewer," % 6F %12G\n",sigma,error);
263: }
264: PetscViewerASCIIPrintf(viewer,"\n");
265: }
266: }
267: return(0);
268: }
272: /*@C
273: SVDCreate - Creates the default SVD context.
275: Collective on MPI_Comm
277: Input Parameter:
278: . comm - MPI communicator
280: Output Parameter:
281: . svd - location to put the SVD context
283: Note:
284: The default SVD type is SVDCROSS
286: Level: beginner
288: .seealso: SVDSetUp(), SVDSolve(), SVDDestroy(), SVD
289: @*/
290: PetscErrorCode SVDCreate(MPI_Comm comm,SVD *outsvd)
291: {
293: SVD svd;
297: *outsvd = 0;
298: SlepcHeaderCreate(svd,_p_SVD,struct _SVDOps,SVD_CLASSID,-1,"SVD","Singular Value Decomposition","SVD",comm,SVDDestroy,SVDView);
300: svd->OP = PETSC_NULL;
301: svd->ip = PETSC_NULL;
302: svd->ds = PETSC_NULL;
303: svd->A = PETSC_NULL;
304: svd->AT = PETSC_NULL;
305: svd->transmode = (SVDTransposeMode)PETSC_DECIDE;
306: svd->sigma = PETSC_NULL;
307: svd->perm = PETSC_NULL;
308: svd->U = PETSC_NULL;
309: svd->V = PETSC_NULL;
310: svd->IS = PETSC_NULL;
311: svd->tl = PETSC_NULL;
312: svd->tr = PETSC_NULL;
313: svd->rand = PETSC_NULL;
314: svd->which = SVD_LARGEST;
315: svd->n = 0;
316: svd->nconv = 0;
317: svd->nsv = 1;
318: svd->ncv = 0;
319: svd->mpd = 0;
320: svd->nini = 0;
321: svd->its = 0;
322: svd->max_it = 0;
323: svd->tol = PETSC_DEFAULT;
324: svd->errest = PETSC_NULL;
325: svd->data = PETSC_NULL;
326: svd->setupcalled = 0;
327: svd->reason = SVD_CONVERGED_ITERATING;
328: svd->numbermonitors = 0;
329: svd->matvecs = 0;
330: svd->trackall = PETSC_FALSE;
332: PetscRandomCreate(comm,&svd->rand);
333: PetscRandomSetSeed(svd->rand,0x12345678);
334: PetscLogObjectParent(svd,svd->rand);
335: *outsvd = svd;
336: return(0);
337: }
338:
341: /*@
342: SVDReset - Resets the SVD context to the setupcalled=0 state and removes any
343: allocated objects.
345: Collective on SVD
347: Input Parameter:
348: . svd - singular value solver context obtained from SVDCreate()
350: Level: advanced
352: .seealso: SVDDestroy()
353: @*/
354: PetscErrorCode SVDReset(SVD svd)
355: {
360: if (svd->ops->reset) { (svd->ops->reset)(svd); }
361: if (svd->ip) { IPReset(svd->ip); }
362: if (svd->ds) { DSReset(svd->ds); }
363: MatDestroy(&svd->OP);
364: MatDestroy(&svd->A);
365: MatDestroy(&svd->AT);
366: VecDestroy(&svd->tl);
367: VecDestroy(&svd->tr);
368: if (svd->n) {
369: PetscFree(svd->sigma);
370: PetscFree(svd->perm);
371: PetscFree(svd->errest);
372: VecDestroyVecs(svd->n,&svd->U);
373: VecDestroyVecs(svd->n,&svd->V);
374: }
375: svd->transmode = (SVDTransposeMode)PETSC_DECIDE;
376: svd->matvecs = 0;
377: svd->setupcalled = 0;
378: return(0);
379: }
383: /*@C
384: SVDDestroy - Destroys the SVD context.
386: Collective on SVD
388: Input Parameter:
389: . svd - singular value solver context obtained from SVDCreate()
391: Level: beginner
393: .seealso: SVDCreate(), SVDSetUp(), SVDSolve()
394: @*/
395: PetscErrorCode SVDDestroy(SVD *svd)
396: {
398:
400: if (!*svd) return(0);
402: if (--((PetscObject)(*svd))->refct > 0) { *svd = 0; return(0); }
403: SVDReset(*svd);
404: PetscObjectDepublish(*svd);
405: if ((*svd)->ops->destroy) { (*(*svd)->ops->destroy)(*svd); }
406: IPDestroy(&(*svd)->ip);
407: DSDestroy(&(*svd)->ds);
408: PetscRandomDestroy(&(*svd)->rand);
409: SVDMonitorCancel(*svd);
410: PetscHeaderDestroy(svd);
411: return(0);
412: }
416: /*@C
417: SVDSetType - Selects the particular solver to be used in the SVD object.
419: Logically Collective on SVD
421: Input Parameters:
422: + svd - the singular value solver context
423: - type - a known method
425: Options Database Key:
426: . -svd_type <method> - Sets the method; use -help for a list
427: of available methods
428:
429: Notes:
430: See "slepc/include/slepcsvd.h" for available methods. The default
431: is SVDCROSS.
433: Normally, it is best to use the SVDSetFromOptions() command and
434: then set the SVD type from the options database rather than by using
435: this routine. Using the options database provides the user with
436: maximum flexibility in evaluating the different available methods.
437: The SVDSetType() routine is provided for those situations where it
438: is necessary to set the iterative solver independently of the command
439: line or options database.
441: Level: intermediate
443: .seealso: SVDType
444: @*/
445: PetscErrorCode SVDSetType(SVD svd,const SVDType type)
446: {
447: PetscErrorCode ierr,(*r)(SVD);
448: PetscBool match;
454: PetscObjectTypeCompare((PetscObject)svd,type,&match);
455: if (match) return(0);
457: PetscFListFind(SVDList,((PetscObject)svd)->comm,type,PETSC_TRUE,(void (**)(void)) &r);
458: if (!r) SETERRQ1(((PetscObject)svd)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown SVD type given: %s",type);
460: if (svd->ops->destroy) { (*svd->ops->destroy)(svd); }
461: PetscMemzero(svd->ops,sizeof(struct _SVDOps));
463: svd->setupcalled = 0;
464: PetscObjectChangeTypeName((PetscObject)svd,type);
465: (*r)(svd);
466: return(0);
467: }
471: /*@C
472: SVDGetType - Gets the SVD type as a string from the SVD object.
474: Not Collective
476: Input Parameter:
477: . svd - the singular value solver context
479: Output Parameter:
480: . name - name of SVD method
482: Level: intermediate
484: .seealso: SVDSetType()
485: @*/
486: PetscErrorCode SVDGetType(SVD svd,const SVDType *type)
487: {
491: *type = ((PetscObject)svd)->type_name;
492: return(0);
493: }
497: /*@C
498: SVDRegister - See SVDRegisterDynamic()
500: Level: advanced
501: @*/
502: PetscErrorCode SVDRegister(const char *sname,const char *path,const char *name,PetscErrorCode (*function)(SVD))
503: {
505: char fullname[PETSC_MAX_PATH_LEN];
508: PetscFListConcat(path,name,fullname);
509: PetscFListAdd(&SVDList,sname,fullname,(void (*)(void))function);
510: return(0);
511: }
515: /*@
516: SVDRegisterDestroy - Frees the list of SVD methods that were
517: registered by SVDRegisterDynamic().
519: Not Collective
521: Level: advanced
523: .seealso: SVDRegisterDynamic(), SVDRegisterAll()
524: @*/
525: PetscErrorCode SVDRegisterDestroy(void)
526: {
530: PetscFListDestroy(&SVDList);
531: SVDRegisterAllCalled = PETSC_FALSE;
532: return(0);
533: }
537: /*@
538: SVDSetIP - Associates an inner product object to the
539: singular value solver.
541: Collective on SVD
543: Input Parameters:
544: + svd - singular value solver context obtained from SVDCreate()
545: - ip - the inner product object
547: Note:
548: Use SVDGetIP() to retrieve the inner product context (for example,
549: to free it at the end of the computations).
551: Level: advanced
553: .seealso: SVDGetIP()
554: @*/
555: PetscErrorCode SVDSetIP(SVD svd,IP ip)
556: {
563: PetscObjectReference((PetscObject)ip);
564: IPDestroy(&svd->ip);
565: svd->ip = ip;
566: PetscLogObjectParent(svd,svd->ip);
567: return(0);
568: }
572: /*@C
573: SVDGetIP - Obtain the inner product object associated
574: to the singular value solver object.
576: Not Collective
578: Input Parameters:
579: . svd - singular value solver context obtained from SVDCreate()
581: Output Parameter:
582: . ip - inner product context
584: Level: advanced
586: .seealso: SVDSetIP()
587: @*/
588: PetscErrorCode SVDGetIP(SVD svd,IP *ip)
589: {
595: if (!svd->ip) {
596: IPCreate(((PetscObject)svd)->comm,&svd->ip);
597: PetscLogObjectParent(svd,svd->ip);
598: }
599: *ip = svd->ip;
600: return(0);
601: }
605: /*@
606: SVDSetDS - Associates a direct solver object to the singular value solver.
608: Collective on SVD
610: Input Parameters:
611: + svd - singular value solver context obtained from SVDCreate()
612: - ds - the direct solver object
614: Note:
615: Use SVDGetDS() to retrieve the direct solver context (for example,
616: to free it at the end of the computations).
618: Level: advanced
620: .seealso: SVDGetDS()
621: @*/
622: PetscErrorCode SVDSetDS(SVD svd,DS ds)
623: {
630: PetscObjectReference((PetscObject)ds);
631: DSDestroy(&svd->ds);
632: svd->ds = ds;
633: PetscLogObjectParent(svd,svd->ds);
634: return(0);
635: }
639: /*@C
640: SVDGetDS - Obtain the direct solver object associated to the singular value
641: solver object.
643: Not Collective
645: Input Parameters:
646: . svd - singular value solver context obtained from SVDCreate()
648: Output Parameter:
649: . ds - direct solver context
651: Level: advanced
653: .seealso: SVDSetDS()
654: @*/
655: PetscErrorCode SVDGetDS(SVD svd,DS *ds)
656: {
662: if (!svd->ds) {
663: DSCreate(((PetscObject)svd)->comm,&svd->ds);
664: PetscLogObjectParent(svd,svd->ds);
665: }
666: *ds = svd->ds;
667: return(0);
668: }